read_surface.f90 Source File


This file depends on

sourcefile~~read_surface.f90~~EfferentGraph sourcefile~read_surface.f90 read_surface.f90 sourcefile~mod_grid_phy_lmdz.f90 mod_grid_phy_lmdz.f90 sourcefile~read_surface.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_para.f90 mod_phys_lmdz_para.f90 sourcefile~read_surface.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~iophy.f90 iophy.F90 sourcefile~read_surface.f90->sourcefile~iophy.f90 sourcefile~dimphy.f90 dimphy.f90 sourcefile~read_surface.f90->sourcefile~dimphy.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_mpi_data.f90 mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_transfert_para.f90 mod_phys_lmdz_transfert_para.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~print_control_mod.f90 print_control_mod.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~print_control_mod.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~iophy.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~iophy.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~iophy.f90->sourcefile~dimphy.f90 sourcefile~phys_output_var_mod.f90 phys_output_var_mod.f90 sourcefile~iophy.f90->sourcefile~phys_output_var_mod.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90 mod_phys_lmdz_omp_transfert.f90 sourcefile~iophy.f90->sourcefile~mod_phys_lmdz_omp_transfert.f90 sourcefile~lmdz_xios.f90 lmdz_xios.F90 sourcefile~iophy.f90->sourcefile~lmdz_xios.f90 sourcefile~wxios_mod.f90 wxios_mod.F90 sourcefile~iophy.f90->sourcefile~wxios_mod.f90 sourcefile~iophy.f90->sourcefile~print_control_mod.f90 sourcefile~clesphys_mod_h.f90 clesphys_mod_h.f90 sourcefile~iophy.f90->sourcefile~clesphys_mod_h.f90 sourcefile~aero_mod.f90 aero_mod.f90 sourcefile~iophy.f90->sourcefile~aero_mod.f90 sourcefile~phys_output_var_mod.f90->sourcefile~dimphy.f90 sourcefile~phys_output_var_mod.f90->sourcefile~clesphys_mod_h.f90 sourcefile~strings_mod.f90 strings_mod.f90 sourcefile~phys_output_var_mod.f90->sourcefile~strings_mod.f90 sourcefile~config_ocean_skin_m.f90 config_ocean_skin_m.F90 sourcefile~phys_output_var_mod.f90->sourcefile~config_ocean_skin_m.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~wxios_mod.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~wxios_mod.f90->sourcefile~dimphy.f90 sourcefile~wxios_mod.f90->sourcefile~lmdz_xios.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~wxios_mod.f90->sourcefile~print_control_mod.f90 sourcefile~geometry_mod.f90 geometry_mod.f90 sourcefile~wxios_mod.f90->sourcefile~geometry_mod.f90 sourcefile~iniprint_mod_h.f90 iniprint_mod_h.f90 sourcefile~wxios_mod.f90->sourcefile~iniprint_mod_h.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~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_mpi_data.f90->sourcefile~print_control_mod.f90 sourcefile~lmdz_mpi.f90 lmdz_mpi.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_mpi.f90 sourcefile~lmdz_cppkeys_wrapper.f90 lmdz_cppkeys_wrapper.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_data.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~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 sourcefile~geometry_mod.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~geometry_mod.f90->sourcefile~nrtype.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~lmdz_mpi.f90 sourcefile~infotrac_phy.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~infotrac_phy.f90->sourcefile~iniprint_mod_h.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~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~strings_mod.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

Contents

Source Code


Source Code

       subroutine read_surface(name,surfa)

     
! common
! ------
       USE ioipsl
!       USE comgeomphy
       USE dimphy
       USE mod_phys_lmdz_para
       USE iophy
       USE netcdf, ONLY: nf90_inq_varid,nf90_noerr,nf90_get_var,nf90_nowrite,nf90_inq_varid,nf90_open
       USE mod_grid_phy_lmdz, ONLY: nbp_lon,  nbp_lat, klon_glo

!!USE paramet_mod_h
IMPLICIT NONE



       character*10 name
       character*10 varname
!
       real tmp_dyn(nbp_lon+1,nbp_lat)
       real tmp_dyn_glo(nbp_lon+1,nbp_lat)
!       real tmp_dyn_glo(nbp_lon,nbp_lat)
       REAL tmp_dyn_invers(nbp_lon+1,nbp_lat)
       real tmp_dyn_invers_glo(nbp_lon+1,nbp_lat)
!       real tmp_dyn_invers_glo(nbp_lon,nbp_lat)
       real tmp_fi(klon)
       real tmp_fi_glo(klon_glo)
       real surfa(klon,5)
       real surfa_glo(klon_glo,5)
!
       integer ncid, varid, rcode, varlatid,tmpid
       integer start_(2),count_(2)
       integer i,j,l,ig
       character*1 str1

!JE20140526<<
      character*4 ::  latstr,aux4s
      logical :: outcycle, isinversed
      real, dimension(nbp_lat) :: lats
      real, dimension(nbp_lat) :: lats_glo
      integer, dimension(1) :: start_j,endj
!JE20140526>>
!$OMP MASTER
       IF (is_mpi_root .AND. is_omp_root) THEN

       print*,'Lecture du fichier donnees_lisa.nc' 
       rcode=nf90_open('donnees_lisa.nc',nf90_nowrite,ncid)
       if ( rcode /= 0 ) then ; call abort_physic('LMDZ','open donnees_lisa.nc dans read_vent',1) ; endif


!JE20140526<<: check if are inversed or not the latitude grid in donnes_lisa
      outcycle=.false.
      latstr='null'
      isinversed=.false.
      do i=1,5
          if (i==1) aux4s='latu'
          if (i==2) aux4s='LATU'
          if (i==3) aux4s='LatU'
          if (i==4) aux4s='Latu'
          if (i==5) aux4s='latU'
          rcode = nf90_inq_varid(ncid, aux4s, tmpid)
          IF ((rcode==0).and.( .not. outcycle )) THEN
            outcycle=.true.
            varlatid=tmpid
          ENDIF
      enddo ! check if it inversed lat
      start_j(1)=1
      endj(1)=nbp_lat
      rcode = nf90_get_var(ncid, varlatid, lats_glo, start_j, endj)
      if ( .not. outcycle ) then ; call abort_physic('LMDZ','get lat dans read_surface',1) ; endif



! check if netcdf is latitude inversed or not. 
      if (lats_glo(1)<lats_glo(2)) isinversed=.true.
! JE20140526>>


       DO i=1,5
          write(str1,'(i1)') i
          varname=trim(name)//str1
          rcode=nf90_inq_varid(ncid,trim(varname),varid)
          if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get'//varname//'  dans read_vent',1) ; endif
!          varid=nf90_inq_varid(ncid,varname,rcode)

!  dimensions pour les champs scalaires et le vent zonal
!  -----------------------------------------------------

          start_(1)=1
          start_(2)=1      
          count_(1)=nbp_lon+1
!          count_(1)=iip1
          count_(2)=nbp_lat
!          count_(2)=jjp1

! mise a zero des tableaux 
! ------------------------
          tmp_dyn(:,:)=0.0
          tmp_fi(:)=0.0
! Lecture
! -----------------------
          rcode = nf90_get_var(ncid, varid, tmp_dyn_glo, start_, count_)
          if ( rcode /= 0 ) then ; call abort_physic('LMDZ','get'//varname//'  dans read_vent',1) ; endif

!      call dump2d(iip1,jjp1,tmp_dyn,'tmp_dyn   ')
       DO j=1, nbp_lat
          DO ig=1, nbp_lon+1
             tmp_dyn_invers_glo(ig,j)=tmp_dyn_glo(ig,nbp_lat-j+1)
          ENDDO
       ENDDO

       
           if (isinversed) then
                        call local_gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
     & tmp_dyn_invers_glo, tmp_fi_glo)
           else      
                        call local_gr_dyn_fi(1, nbp_lon+1, nbp_lat, klon_glo, &
     &   tmp_dyn_glo, tmp_fi_glo)
           endif
!JE20140526>>
!
          DO j=1,klon_glo

                surfa_glo(j,i)=tmp_fi_glo(j)

          ENDDO ! Fin de recopie du tableau
!
       ENDDO ! Fin boucle 1 a 5
       print*,'Passage Grille Dyn -> Phys'


      ENDIF !mpi 
!$OMP END MASTER
!$OMP BARRIER
      call scatter(surfa_glo,surfa)


       return
       end subroutine read_surface

      SUBROUTINE local_gr_dyn_fi(nfield,im,jm,ngrid,pdyn,pfi)
      IMPLICIT NONE
!=======================================================================
!   passage d'un champ de la grille scalaire a la grille physique
!=======================================================================

!-----------------------------------------------------------------------
!   declarations:
!   -------------

      INTEGER im,jm,ngrid,nfield
      REAL pdyn(im,jm,nfield)
      REAL pfi(ngrid,nfield)

      INTEGER j,ifield,ig

!-----------------------------------------------------------------------
!   calcul:
!   -------

      IF(ngrid.NE.2+(jm-2)*(im-1).AND.ngrid.NE.1)                          &
     &    STOP 'probleme de dim'
!   traitement des poles
      CALL SCOPY(nfield,pdyn,im*jm,pfi,ngrid)
      CALL SCOPY(nfield,pdyn(1,jm,1),im*jm,pfi(ngrid,1),ngrid)

!   traitement des point normaux
      DO ifield=1,nfield
         DO j=2,jm-1
            ig=2+(j-2)*(im-1)
            CALL SCOPY(im-1,pdyn(1,j,ifield),1,pfi(ig,ifield),1)
         ENDDO
      ENDDO

      RETURN
      END SUBROUTINE local_gr_dyn_fi