divgrad2_loc.f90 Source File


This file depends on

sourcefile~~divgrad2_loc.f90~~EfferentGraph sourcefile~divgrad2_loc.f90 divgrad2_loc.f90 sourcefile~parallel_lmdz.f90 parallel_lmdz.F90 sourcefile~divgrad2_loc.f90->sourcefile~parallel_lmdz.f90 sourcefile~paramet_mod_h.f90 paramet_mod_h.f90 sourcefile~divgrad2_loc.f90->sourcefile~paramet_mod_h.f90 sourcefile~comgeom2_mod_h.f90 comgeom2_mod_h.f90 sourcefile~divgrad2_loc.f90->sourcefile~comgeom2_mod_h.f90 sourcefile~divgrad2_mod.f90 divgrad2_mod.f90 sourcefile~divgrad2_loc.f90->sourcefile~divgrad2_mod.f90 sourcefile~mod_hallo.f90 mod_hallo.f90 sourcefile~divgrad2_loc.f90->sourcefile~mod_hallo.f90 sourcefile~times.f90 times.f90 sourcefile~divgrad2_loc.f90->sourcefile~times.f90 sourcefile~comdissipn_mod_h.f90 comdissipn_mod_h.f90 sourcefile~divgrad2_loc.f90->sourcefile~comdissipn_mod_h.f90 sourcefile~parallel_lmdz.f90->sourcefile~paramet_mod_h.f90 sourcefile~vampir.f90 vampir.F90 sourcefile~parallel_lmdz.f90->sourcefile~vampir.f90 sourcefile~lmdz_mpi.f90 lmdz_mpi.F90 sourcefile~parallel_lmdz.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_const_mpi.f90 mod_const_mpi.f90 sourcefile~parallel_lmdz.f90->sourcefile~mod_const_mpi.f90 sourcefile~iniprint_mod_h.f90 iniprint_mod_h.f90 sourcefile~parallel_lmdz.f90->sourcefile~iniprint_mod_h.f90 sourcefile~control_mod.f90 control_mod.f90 sourcefile~parallel_lmdz.f90->sourcefile~control_mod.f90 sourcefile~wxios_mod.f90 wxios_mod.F90 sourcefile~parallel_lmdz.f90->sourcefile~wxios_mod.f90 sourcefile~comgeom2_mod_h.f90->sourcefile~paramet_mod_h.f90 sourcefile~divgrad2_mod.f90->sourcefile~parallel_lmdz.f90 sourcefile~allocate_field_mod.f90 allocate_field_mod.f90 sourcefile~divgrad2_mod.f90->sourcefile~allocate_field_mod.f90 sourcefile~bands.f90 bands.f90 sourcefile~divgrad2_mod.f90->sourcefile~bands.f90 sourcefile~mod_hallo.f90->sourcefile~parallel_lmdz.f90 sourcefile~mod_hallo.f90->sourcefile~paramet_mod_h.f90 sourcefile~mod_hallo.f90->sourcefile~lmdz_mpi.f90 sourcefile~times.f90->sourcefile~parallel_lmdz.f90 sourcefile~times.f90->sourcefile~paramet_mod_h.f90 sourcefile~times.f90->sourcefile~lmdz_mpi.f90 sourcefile~allocate_field_mod.f90->sourcefile~parallel_lmdz.f90 sourcefile~allocate_field_mod.f90->sourcefile~paramet_mod_h.f90 sourcefile~allocate_field_mod.f90->sourcefile~mod_hallo.f90 sourcefile~bands.f90->sourcefile~parallel_lmdz.f90 sourcefile~bands.f90->sourcefile~times.f90 sourcefile~mod_phys_lmdz_para.f90 mod_phys_lmdz_para.f90 sourcefile~bands.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~lmdz_cppkeys_wrapper.f90 lmdz_cppkeys_wrapper.F90 sourcefile~bands.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~wxios_mod.f90->sourcefile~iniprint_mod_h.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~dimphy.f90 dimphy.f90 sourcefile~wxios_mod.f90->sourcefile~dimphy.f90 sourcefile~lmdz_xios.f90 lmdz_xios.F90 sourcefile~wxios_mod.f90->sourcefile~lmdz_xios.f90 sourcefile~mod_phys_lmdz_mpi_data.f90 mod_phys_lmdz_mpi_data.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~geometry_mod.f90 geometry_mod.f90 sourcefile~wxios_mod.f90->sourcefile~geometry_mod.f90 sourcefile~infotrac_phy.f90 infotrac_phy.F90 sourcefile~wxios_mod.f90->sourcefile~infotrac_phy.f90 sourcefile~strings_mod.f90 strings_mod.f90 sourcefile~wxios_mod.f90->sourcefile~strings_mod.f90 sourcefile~print_control_mod.f90 print_control_mod.f90 sourcefile~wxios_mod.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_transfert_para.f90 mod_phys_lmdz_transfert_para.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~mod_grid_phy_lmdz.f90 mod_grid_phy_lmdz.f90 sourcefile~wxios_mod.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~nrtype.f90 nrtype.f90 sourcefile~wxios_mod.f90->sourcefile~nrtype.f90 sourcefile~ioipsl_getin_p_mod.f90 ioipsl_getin_p_mod.f90 sourcefile~wxios_mod.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_omp_data.f90 mod_phys_lmdz_omp_data.F90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~geometry_mod.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~geometry_mod.f90->sourcefile~nrtype.f90 sourcefile~infotrac_phy.f90->sourcefile~iniprint_mod_h.f90 sourcefile~infotrac_phy.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~infotrac_phy.f90->sourcefile~strings_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~readtracfiles_mod.f90 readTracFiles_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~readtracfiles_mod.f90 sourcefile~lmdz_reprobus_wrappers.f90 lmdz_reprobus_wrappers.F90 sourcefile~infotrac_phy.f90->sourcefile~lmdz_reprobus_wrappers.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90 mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90 mod_phys_lmdz_mpi_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_transfert.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~strings_mod.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~readtracfiles_mod.f90->sourcefile~strings_mod.f90 sourcefile~readtracfiles_mod.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~lmdz_reprobus_wrappers.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~dimphy.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~print_control_mod.f90

Contents

Source Code


Source Code

SUBROUTINE divgrad2_loc ( klevel, h, deltapres, lh, divgra_out )
  !
  ! P. Le Van
  !
  !   ***************************************************************
  !
  ! .....   calcul de  (div( grad ))   de (  pext * h ) .....
  !   ****************************************************************
  !   h ,klevel,lh et pext  sont des arguments  d'entree pour le s-prg
  !     divgra     est  un argument  de sortie pour le s-prg
  !
  USE comgeom2_mod_h
  USE comdissipn_mod_h
  USE parallel_lmdz
  USE times
  USE mod_hallo
  USE divgrad2_mod
  USE dimensions_mod, ONLY: iim, jjm, llm, ndm
USE paramet_mod_h
IMPLICIT NONE
  !



  !    .......    variables en arguments   .......
  !
  INTEGER :: klevel
  REAL :: h( ijb_u:ije_u,klevel ), deltapres( ijb_u:ije_u,klevel )
  REAL :: divgra_out( ijb_u:ije_u,klevel)
  !    .......    variables  locales    ..........
  !
  REAL :: signe, nudivgrs, sqrtps( ijb_u:ije_u,llm )
  INTEGER :: l,ij,iter,lh
  !    ...................................................................
  Type(Request),SAVE :: request_dissip
!$OMP THREADPRIVATE(request_dissip)
  INTEGER :: ijb,ije

  !
  !
  signe    = (-1.)**lh
  nudivgrs = signe * cdivh

   ! CALL SCOPY ( ip1jmp1 * klevel, h, 1, divgra, 1 )
  ijb=ij_begin
  ije=ij_end
!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO l = 1, klevel
    divgra(ijb:ije,l)=h(ijb:ije,l)
  ENDDO
!$OMP END DO NOWAIT
  !
!$OMP BARRIER
   call Register_Hallo_u(divgra,llm,1,1,1,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)
!$OMP BARRIER

  CALL laplacien_loc( klevel, divgra, divgra )

!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO l = 1, klevel
   DO ij = ijb, ije
    sqrtps( ij,l ) = SQRT( deltapres(ij,l) )
   ENDDO
  ENDDO
!$OMP END DO NOWAIT

  !
!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO l = 1, klevel
    DO ij = ijb, ije
     divgra(ij,l) = divgra(ij,l) * sqrtps(ij,l)
    ENDDO
  ENDDO
!$OMP END DO NOWAIT

  !    ........    Iteration de l'operateur  laplacien_gam    ........
  !
  DO  iter = 1, lh - 2
!$OMP BARRIER
   call Register_Hallo_u(divgra,llm,1,1,1,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)

!$OMP BARRIER


   CALL laplacien_gam_loc(klevel,cuvscvgam2,cvuscugam2,unsair_gam2, &
         unsapolnga2, unsapolsga2,  divgra, divgra )
  ENDDO
  !
  !    ...............................................................

!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO l = 1, klevel
    DO ij = ijb, ije
      divgra(ij,l) = divgra(ij,l) * sqrtps(ij,l)
    ENDDO
  ENDDO
!$OMP END DO NOWAIT
  !
!$OMP BARRIER
   call Register_Hallo_u(divgra,llm,1,1,1,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)
!$OMP BARRIER

  CALL laplacien_loc ( klevel, divgra, divgra )
  !
!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO l  = 1,klevel
  DO ij = ijb,ije
  divgra_out(ij,l) =  nudivgrs * divgra(ij,l) / deltapres(ij,l)
  ENDDO
  ENDDO
!$OMP END DO NOWAIT

  RETURN
END SUBROUTINE divgrad2_loc