tab_Jf_tcheby Subroutine

public subroutine tab_Jf_tcheby(nx1, nx2, nb_pts, nvarx, nvary, nb_var, tab_tche1, tab_tche2, tab_Jf, imask, multi_thread)

Tableau des dérivées par rapport aux coefficients de tab_tche

Arguments

Type IntentOptional Attributes Name
integer(kind=I4), intent(in) :: nx1

nbre de points de calcul selon x

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

nbre de points de calcul selon y

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

nbre de points

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

degré max de Tchebychev utilisé selon x

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

degré max de Tchebychev utilisé selon y

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

nombre de fonctions de base utilisées

real(kind=R8), intent(in), dimension(1:nx1, 1:nvarx+1) :: tab_tche1

tableau des valeurs calculées

real(kind=R8), intent(in), dimension(1:nx2, 1:nvary+1) :: tab_tche2

tableau des valeurs calculées

real(kind=R8), intent(out), allocatable, dimension(:,:) :: tab_Jf

tableau des dérivées

integer(kind=I4), intent(in), optional, dimension(1:nx1, 1:nx2) :: imask

masque

logical(kind=I4), intent(in), optional :: multi_thread

Calls

proc~~tab_jf_tcheby~~CallsGraph proc~tab_jf_tcheby tab_Jf_tcheby omp_get_num_procs omp_get_num_procs proc~tab_jf_tcheby->omp_get_num_procs

Called by

proc~~tab_jf_tcheby~~CalledByGraph proc~tab_jf_tcheby tab_Jf_tcheby proc~least_squares_tcheby least_squares_tcheby proc~least_squares_tcheby->proc~tab_jf_tcheby program~test_tchebychev test_tchebychev program~test_tchebychev->proc~least_squares_tcheby

Source Code

   subroutine tab_Jf_tcheby(nx1, nx2, nb_pts, nvarx, nvary, nb_var, tab_tche1, tab_tche2, tab_Jf, imask, multi_thread)
   !! Tableau des dérivées par rapport aux coefficients de tab_tche
   implicit none
   integer(kind=I4), intent(in ) :: nb_var                                    !! *nombre de fonctions de base utilisées*
   integer(kind=I4), intent(in ) :: nx1                                       !! *nbre de points de calcul selon x*
   integer(kind=I4), intent(in ) :: nx2                                       !! *nbre de points de calcul selon y*
   integer(kind=I4), intent(in ) :: nb_pts                                    !! *nbre de points*
   integer(kind=I4), intent(in ) :: nvarx                                     !! *degré max de Tchebychev utilisé selon x*
   integer(kind=I4), intent(in ) :: nvary                                     !! *degré max de Tchebychev utilisé selon y*
   real(kind=R8),    intent(in ), dimension(1:nx1, 1:nvarx+1) :: tab_tche1    !! *tableau des valeurs calculées*
   real(kind=R8),    intent(in ), dimension(1:nx2, 1:nvary+1) :: tab_tche2    !! *tableau des valeurs calculées*
   integer(kind=I4), intent(in ), dimension(1:nx1, 1:nx2), optional :: imask  !! *masque*
   real(kind=R8),    intent(out), allocatable, dimension(:,:) :: tab_Jf       !! *tableau des dérivées*
   logical(kind=I4), intent(in ), optional :: multi_thread

      integer(kind=I4) :: ivar, jvar, ij, ipt, jpt, ijpt, ibatch, nb_threads
      logical(kind=I4) :: lmask
      logical(kind=I4) :: mlth

      mlth = .false.
      if ( present( multi_thread ) ) mlth = multi_thread

      nb_threads = omp_get_num_procs()
      ibatch = max( nx2/nb_threads, 1 )

      lmask = present(imask)

      allocate( tab_Jf(1:nb_pts, 1:nb_var) )
      ! table des dérivées des poly par rapport aux coeff, c'est donc UN produit de polynômes

      !$OMP PARALLEL DEFAULT(SHARED) NUM_THREADS(nb_threads) IF(mlth)
      !$OMP DO SCHEDULE (STATIC,ibatch) PRIVATE(ijpt, ipt, jpt, ij, ivar, jvar)

      do jpt = 1, nx2
         do ipt = 1, nx1

            if ( lmask ) then
               if ( imask(ipt, jpt)==0 ) cycle
            endif

            ijpt = (jpt-1)*nx1 + ipt

            ij = 0
            do jvar = 0, nvary
               do ivar = 0, nvarx
                  ij = ij +1
                  tab_Jf(ijpt, ij) = tab_tche1(ipt, ivar +1)*tab_tche2(jpt, jvar +1)
               enddo
            enddo
         enddo
      enddo

      !$OMP END DO
      !$OMP END PARALLEL

   return
   endsubroutine tab_Jf_tcheby