coeff_tcheby_vers_monome Subroutine

public subroutine coeff_tcheby_vers_monome(coeff_t, coeff_m, deg)

Note

Transformation des coefficients de Tchebychev en coefficients de monômes

Un monôme de degré p reçoit de Ti(x) (polynôme de Tcheby de degré i) :

Arguments

Type IntentOptional Attributes Name
real(kind=R8), intent(in), dimension(0:deg) :: coeff_t

coefficients de la CL de polynômes de Tchebychev

real(kind=R8), intent(out), dimension(0:deg) :: coeff_m

coefficients de la CL de monômes

integer(kind=I4), intent(in) :: deg

degré du polynôme


Called by

proc~~coeff_tcheby_vers_monome~~CalledByGraph proc~coeff_tcheby_vers_monome coeff_tcheby_vers_monome proc~coeff_tcheby_xy_vers_monome coeff_tcheby_xy_vers_monome proc~coeff_tcheby_xy_vers_monome->proc~coeff_tcheby_vers_monome proc~genere_surf_poly genere_surf_poly proc~genere_surf_poly->proc~coeff_tcheby_vers_monome proc~coeff_poly_tcheby_xy_vers_poly_monome coeff_poly_tcheby_xy_vers_poly_monome proc~coeff_poly_tcheby_xy_vers_poly_monome->proc~coeff_tcheby_xy_vers_monome program~test_tchebychev test_tchebychev program~test_tchebychev->proc~genere_surf_poly proc~least_squares_tcheby least_squares_tcheby program~test_tchebychev->proc~least_squares_tcheby proc~least_squares_tcheby->proc~coeff_poly_tcheby_xy_vers_poly_monome

Source Code

   subroutine coeff_tcheby_vers_monome(coeff_t, coeff_m, deg)
   implicit none
   integer(kind=I4), intent(in) :: deg                      !! *degré du polynôme*
   real(kind=R8), intent(in ), dimension(0:deg) :: coeff_t  !! *coefficients de la CL de polynômes de Tchebychev*
   real(kind=R8), intent(out), dimension(0:deg) :: coeff_m  !! *coefficients de la CL de monômes*

      integer(kind=I4) :: k, n, p, q
      real(kind=R8) :: tkm2q, tmp

      n = deg
      coeff_m(0:deg) = 0._R8

      select case (n)

         case (0)
            coeff_m(0)   =   coeff_t(0)
         case (1)
            coeff_m(0:1) = [ coeff_t(0)             , coeff_t(1) ]
         case (2)
            coeff_m(0:2) = [ coeff_t(0) - coeff_t(2), coeff_t(1)                 , 2 * coeff_t(2) ]
         case (3)
            coeff_m(0:3) = [ coeff_t(0) - coeff_t(2), coeff_t(1) - 3 * coeff_t(3), 2 * coeff_t(2), 4 * coeff_t(3) ]
         case (4:)
            coeff_m(0:3) = [ coeff_t(0) - coeff_t(2), coeff_t(1) - 3 * coeff_t(3), 2 * coeff_t(2), 4 * coeff_t(3) ]

            do k = 4, n

               q = 0
               coeff_m(k-2*q) = coeff_m(k-2*q) + coeff_t(k) * (2**(k-1))

               q = 1
               coeff_m(k-2*q) = coeff_m(k-2*q) + coeff_t(k) * (-k/2._R8) * (2**(k-2))

               do q = 2, k/2

                  tmp = 1._R8
                  do p = 1, q - 1
                     tmp = tmp * ( (k - q - p) / (1._R8*p) )
                  enddo

                  tkm2q = (k/2._R8) * ((-1)**q) * (2**(k-2*q)) * (1._R8/q) * tmp

                  coeff_m(k-2*q) = coeff_m(k-2*q) + coeff_t(k) * tkm2q

               enddo

            enddo
      endselect

   return
   endsubroutine coeff_tcheby_vers_monome