gradiv2_loc.f90 Source File


This file depends on

sourcefile~~gradiv2_loc.f90~~EfferentGraph sourcefile~gradiv2_loc.f90 gradiv2_loc.f90 sourcefile~parallel_lmdz.f90 parallel_lmdz.F90 sourcefile~gradiv2_loc.f90->sourcefile~parallel_lmdz.f90 sourcefile~gradiv2_mod.f90 gradiv2_mod.f90 sourcefile~gradiv2_loc.f90->sourcefile~gradiv2_mod.f90 sourcefile~paramet_mod_h.f90 paramet_mod_h.f90 sourcefile~gradiv2_loc.f90->sourcefile~paramet_mod_h.f90 sourcefile~write_field_p.f90 write_field_p.f90 sourcefile~gradiv2_loc.f90->sourcefile~write_field_p.f90 sourcefile~comgeom_mod_h.f90 comgeom_mod_h.f90 sourcefile~gradiv2_loc.f90->sourcefile~comgeom_mod_h.f90 sourcefile~mod_hallo.f90 mod_hallo.f90 sourcefile~gradiv2_loc.f90->sourcefile~mod_hallo.f90 sourcefile~mod_filtreg_p.f90 mod_filtreg_p.F90 sourcefile~gradiv2_loc.f90->sourcefile~mod_filtreg_p.f90 sourcefile~times.f90 times.f90 sourcefile~gradiv2_loc.f90->sourcefile~times.f90 sourcefile~comdissipn_mod_h.f90 comdissipn_mod_h.f90 sourcefile~gradiv2_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~gradiv2_mod.f90->sourcefile~parallel_lmdz.f90 sourcefile~allocate_field_mod.f90 allocate_field_mod.f90 sourcefile~gradiv2_mod.f90->sourcefile~allocate_field_mod.f90 sourcefile~bands.f90 bands.f90 sourcefile~gradiv2_mod.f90->sourcefile~bands.f90 sourcefile~write_field_p.f90->sourcefile~parallel_lmdz.f90 sourcefile~write_field.f90 write_field.f90 sourcefile~write_field_p.f90->sourcefile~write_field.f90 sourcefile~comgeom_mod_h.f90->sourcefile~paramet_mod_h.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~mod_filtreg_p.f90->sourcefile~parallel_lmdz.f90 sourcefile~mod_filtreg_p.f90->sourcefile~paramet_mod_h.f90 sourcefile~filtreg_mod.f90 filtreg_mod.F90 sourcefile~mod_filtreg_p.f90->sourcefile~filtreg_mod.f90 sourcefile~coefils_mod_h.f90 coefils_mod_h.f90 sourcefile~mod_filtreg_p.f90->sourcefile~coefils_mod_h.f90 sourcefile~mod_filtre_fft_loc.f90 mod_filtre_fft_loc.F90 sourcefile~mod_filtreg_p.f90->sourcefile~mod_filtre_fft_loc.f90 sourcefile~timer_filtre.f90 timer_filtre.f90 sourcefile~mod_filtreg_p.f90->sourcefile~timer_filtre.f90 sourcefile~times.f90->sourcefile~parallel_lmdz.f90 sourcefile~times.f90->sourcefile~paramet_mod_h.f90 sourcefile~times.f90->sourcefile~lmdz_mpi.f90 sourcefile~filtreg_mod.f90->sourcefile~paramet_mod_h.f90 sourcefile~filtreg_mod.f90->sourcefile~comgeom_mod_h.f90 sourcefile~filtreg_mod.f90->sourcefile~coefils_mod_h.f90 sourcefile~comconst_mod.f90 comconst_mod.f90 sourcefile~filtreg_mod.f90->sourcefile~comconst_mod.f90 sourcefile~serre_mod.f90 serre_mod.f90 sourcefile~filtreg_mod.f90->sourcefile~serre_mod.f90 sourcefile~logic_mod.f90 logic_mod.f90 sourcefile~filtreg_mod.f90->sourcefile~logic_mod.f90 sourcefile~mod_fft.f90 mod_fft.F90 sourcefile~mod_filtre_fft_loc.f90->sourcefile~mod_fft.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~strings_mod.f90 strings_mod.f90 sourcefile~write_field.f90->sourcefile~strings_mod.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~lmdz_xios.f90 lmdz_xios.F90 sourcefile~wxios_mod.f90->sourcefile~lmdz_xios.f90 sourcefile~dimphy.f90 dimphy.f90 sourcefile~wxios_mod.f90->sourcefile~dimphy.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~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_fft_wrapper.f90 mod_fft_wrapper.f90 sourcefile~mod_fft.f90->sourcefile~mod_fft_wrapper.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~lmdz_cppkeys_wrapper.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~print_control_mod.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~lmdz_cppkeys_wrapper.f90 sourcefile~infotrac_phy.f90->sourcefile~strings_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~ioipsl_getin_p_mod.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 gradiv2_loc(klevel, xcov, ycov, ld, gdx_out, gdy_out )
  !
  ! P. Le Van
  !
  !   **********************************************************
  !                            ld
  !   calcul  de  (grad (div) )   du vect. v ....
  !
  ! xcov et ycov etant les composant.covariantes de v
  !   **********************************************************
  ! xcont , ycont et ld  sont des arguments  d'entree pour le s-prog
  !  gdx   et  gdy       sont des arguments de sortie pour le s-prog
  !
  !
  USE comgeom_mod_h
  USE comdissipn_mod_h
  USE parallel_lmdz
  USE times
  USE Write_field_p
  USE mod_hallo
  USE mod_filtreg_p
  USE gradiv2_mod
  USE dimensions_mod, ONLY: iim, jjm, llm, ndm
USE paramet_mod_h
IMPLICIT NONE
  !


  !
  ! ........    variables en arguments      ........

  INTEGER :: klevel
  REAL :: xcov( ijb_u:ije_u,klevel ), ycov( ijb_v:ije_v,klevel )
  REAL :: gdx_out( ijb_u:ije_u,klevel ), gdy_out( ijb_v:ije_v,klevel)
  !
  ! ........       variables locales       .........
  !
  REAL      :: tmp_div2(ijb_u:ije_u,llm)
  REAL :: signe, nugrads
  INTEGER :: l,ij,iter,ld
  INTEGER :: ijb,ije,jjb,jje
  Type(Request),SAVE  :: request_dissip
!$OMP THREADPRIVATE(request_dissip)
  !    ........................................................
  !
  !
  !  CALL SCOPY( ip1jmp1 * klevel, xcov, 1, gdx, 1 )
  !  CALL SCOPY(   ip1jm * klevel, ycov, 1, gdy, 1 )

  ijb=ij_begin
  ije=ij_end

!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO   l = 1, klevel
    gdx(ijb:ije,l)=xcov(ijb:ije,l)
  ENDDO
!$OMP END DO NOWAIT

  ijb=ij_begin
  ije=ij_end
  if(pole_sud) ije=ij_end-iip1

!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  DO   l = 1, klevel
    gdy(ijb:ije,l)=ycov(ijb:ije,l)
  ENDDO
!$OMP END DO NOWAIT

!$OMP BARRIER
   call Register_Hallo_v(gdy,llm,1,0,0,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)
!$OMP BARRIER
  !
  !
  signe   = (-1.)**ld
  nugrads = signe * cdivu
  !


  CALL    divergf_loc( klevel, gdx,   gdy , div )
   ! call write_field3d_p('div1',reshape(div,(/iip1,jjp1,llm/)))

  IF( ld.GT.1 )   THEN
!$OMP BARRIER
   call Register_Hallo_u(div,llm,1,1,1,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)
!$OMP BARRIER
    CALL laplacien_loc( klevel, div,  div     )

  !    ......  Iteration de l'operateur laplacien_gam   .......
      ! call write_field3d_p('div2',reshape(div,(/iip1,jjp1,llm/)))

    DO iter = 1, ld -2
!$OMP BARRIER
   call Register_Hallo_u(div,llm,1,1,1,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)

!$OMP BARRIER

     CALL laplacien_gam_loc(klevel,cuvscvgam1,cvuscugam1, &
           unsair_gam1,unsapolnga1, unsapolsga1, &
           div, div       )
    ENDDO
     ! call write_field3d_p('div3',reshape(div,(/iip1,jjp1,llm/)))
  ENDIF

   jjb=jj_begin
   jje=jj_end

   CALL filtreg_p( div   ,jjb_u,jje_u,jjb,jje, jjp1, &
         klevel, 2, 1, .TRUE., 1 )
    ! call exchange_Hallo(div,ip1jmp1,llm,0,1)
!$OMP BARRIER
   call Register_Hallo_u(div,llm,1,1,1,1,Request_dissip)
   call SendRequest(Request_dissip)
!$OMP BARRIER
   call WaitRequest(Request_dissip)

!$OMP BARRIER


   CALL  grad_loc( klevel,  div,   gdx,  gdy )

  !
  ijb=ij_begin
  ije=ij_end

!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
   DO   l = 1, klevel

     if (pole_sud) ije=ij_end
     DO  ij = ijb, ije
      gdx_out( ij,l ) = gdx( ij,l ) * nugrads
     ENDDO

     if (pole_sud) ije=ij_end-iip1
     DO  ij = ijb, ije
      gdy_out( ij,l ) = gdy( ij,l ) * nugrads
     ENDDO

   ENDDO
!$OMP END DO NOWAIT
  !
   RETURN
END SUBROUTINE gradiv2_loc