source: LMDZ6/trunk/libf/dynphy_lonlat/phylmd/ce0l.F90 @ 5408

Last change on this file since 5408 was 5285, checked in by abarral, 7 weeks ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 9.5 KB
RevLine 
[2293]1PROGRAM ce0l
[2336]2!
3!-------------------------------------------------------------------------------
4! Purpose: Initial states and boundary conditions files creation:
5!     * start.nc    for dynamics    (using etat0dyn     routine)
6!     * startphy.nc for physics     (using etat0phys    routine)
7!     * limit.nc    for forced runs (using limit_netcdf routine)
8!-------------------------------------------------------------------------------
9! Notes:
10!     * extrap=.T. (default) for data extrapolation, like for the SSTs when file
11!                   does contain ocean points only.
12!     * "masque" can be:
13!       - read from file "o2a.nc"          (for coupled runs).
[2665]14!       - read from file "startphy0.nc"    (from a previous run).
[2336]15!       - created in etat0phys or etat0dyn (for forced  runs).
16!     It is then passed to limit_netcdf to ensure consistancy.
17!-------------------------------------------------------------------------------
[2293]18  USE ioipsl, ONLY: ioconf_calendar, getin, flininfo, flinopen, flinget, flinclo
[2336]19  USE control_mod,    ONLY: day_step, dayref, nsplit_phys
20  USE etat0dyn,       ONLY: etat0dyn_netcdf
21  USE etat0phys,      ONLY: etat0phys_netcdf
22  USE limit,          ONLY: limit_netcdf
[5084]23  USE netcdf,         ONLY: NF90_OPEN, NF90_NOWRITE, NF90_CLOSE, NF90_NOERR,    &
[2665]24         NF90_INQUIRE_DIMENSION, NF90_INQ_DIMID, NF90_INQ_VARID, NF90_GET_VAR
[4389]25  USE infotrac,       ONLY: init_infotrac
[2336]26  USE dimphy,         ONLY: klon
[2293]27  USE test_disvert_m, ONLY: test_disvert
[2336]28  USE filtreg_mod,    ONLY: inifilr
[2349]29  USE iniphysiq_mod,  ONLY: iniphysiq
[2351]30  USE mod_const_mpi,  ONLY: comm_lmdz
[3815]31
[2336]32#ifdef CPP_PARA
33  USE mod_const_mpi,  ONLY: init_const_mpi
[3815]34  USE parallel_lmdz,  ONLY: init_parallel, mpi_rank, omp_rank, using_mpi
[2336]35  USE bands,          ONLY: read_distrib, distrib_phys
36  USE mod_hallo,      ONLY: init_mod_hallo
37  USE mod_interface_dyn_phys, ONLY: init_interface_dyn_phys
[4689]38  USE lmdz_xios, only: using_xios, xios_finalize
[2336]39#endif
[3815]40
[5282]41  USE iniprint_mod_h
[5281]42  USE comgeom2_mod_h
[2597]43  USE comconst_mod, ONLY: cpp, daysec, dtphys, dtvr, g, kappa, omeg, r, rad, &
44                          pi, jmp1
[2603]45  USE logic_mod, ONLY: iflag_phys, ok_etat0, ok_limit
[2600]46  USE comvert_mod, ONLY: pa, preff, pressure_exner
[2601]47  USE temps_mod, ONLY: calend, day_ini, dt
[4600]48  USE lmdz_mpi
[524]49
[5271]50  USE dimensions_mod, ONLY: iim, jjm, llm, ndm
[5285]51USE paramet_mod_h
[5271]52IMPLICIT NONE
[2293]53
[2336]54!-------------------------------------------------------------------------------
55! Local variables:
[5271]56
[5272]57
[3815]58 
[2336]59  REAL               :: masque(iip1,jjp1)             !--- CONTINENTAL MASK
60  REAL               :: phis  (iip1,jjp1)             !--- GROUND GEOPOTENTIAL
[2293]61  CHARACTER(LEN=256) :: modname, fmt, calnd           !--- CALENDAR TYPE
62  LOGICAL            :: use_filtre_fft
[2336]63  LOGICAL, PARAMETER :: extrap=.FALSE.
[2293]64
[2336]65!--- Local variables for ocean mask reading:
[2293]66  INTEGER            :: nid_o2a, iml_omask, jml_omask, j
67  INTEGER            :: fid, iret, llm_tmp, ttm_tmp, itaul(1)
[2336]68  REAL, ALLOCATABLE  :: lon_omask(:,:), dlon_omask(:), ocemask(:,:)
69  REAL, ALLOCATABLE  :: lat_omask(:,:), dlat_omask(:), ocetmp (:,:)
[2293]70  REAL               :: date, lev(1)
[2665]71
72!--- Local variables for land mask from startphy0 file reading
73  INTEGER            :: nid_sta, nid_nph, nid_msk, nphys
74  REAL, ALLOCATABLE  :: masktmp(:)
75
[3815]76#ifdef CPP_PARA
77  integer ierr
78#else
[2353]79! for iniphysiq in serial mode
80  INTEGER,PARAMETER :: mpi_rank=0
81  INTEGER :: distrib_phys(mpi_rank:mpi_rank)=(jjm-1)*iim+2
82#endif
[2336]83!-------------------------------------------------------------------------------
[2293]84  modname="ce0l"
85
[2336]86!--- Constants
[2293]87  pi     = 4. * ATAN(1.)
88  rad    = 6371229.
89  daysec = 86400.
90  omeg   = 2.*pi/daysec
91  g      = 9.8
92  kappa  = 0.2857143
93  cpp    = 1004.70885
94  jmp1   = jjm + 1
95  preff   = 101325.
96  pa      = 50000.
97
[2221]98  CALL conf_gcm( 99, .TRUE. )
[2331]99  dtvr = daysec/REAL(day_step)
[2336]100  WRITE(lunout,*)'dtvr',dtvr
[2293]101  CALL iniconst()
102  CALL inigeom()
[822]103
[2336]104!--- Calendar choice
[2293]105  calnd='gregorian'
[1319]106  SELECT CASE(calend)
[4361]107    CASE('earth_360d');CALL ioconf_calendar('360_day');   calnd='with 360 days/year'
[2336]108    CASE('earth_365d');CALL ioconf_calendar('noleap'); calnd='with no leap year'
109    CASE('earth_366d');CALL ioconf_calendar('366d');   calnd='with leap years only'
110    CASE('gregorian'); CALL ioconf_calendar('gregorian')
111    CASE('standard');  CALL ioconf_calendar('gregorian')
112    CASE('julian');    CALL ioconf_calendar('julian'); calnd='julian'
113    CASE('proleptic_gregorian'); CALL ioconf_calendar('gregorian')
114  !--- DC Bof...  => IOIPSL a mettre a jour: proleptic_gregorian /= gregorian
115    CASE DEFAULT
116      CALL abort_gcm('ce0l','Bad choice for calendar',1)
[1319]117  END SELECT
[2336]118  WRITE(lunout,*)'CHOSEN CALENDAR: Earth '//TRIM(calnd)
[1279]119
[5267]120
[2336]121#ifdef CPP_PARA
122!--- Physical grid + parallel initializations
123  CALL init_const_mpi()
124  CALL init_parallel()
125  CALL read_distrib()
126  CALL init_mod_hallo()
127#endif
128  WRITE(lunout,*)'---> klon=',klon
129
130!--- Tracers initializations
[4325]131  CALL init_infotrac()
[2336]132
133  CALL inifilr()
[2351]134  CALL iniphysiq(iim,jjm,llm, &
135                 distrib_phys(mpi_rank),comm_lmdz, &
136                 daysec,day_ini,dtphys/nsplit_phys, &
[2349]137                 rlatu,rlatv,rlonu,rlonv,aire,cu,cv,rad,g,r,cpp,iflag_phys)
[2336]138  IF(pressure_exner) CALL test_disvert
139
140#ifdef CPP_PARA
141  IF (mpi_rank==0.AND.omp_rank==0) THEN
142#endif
[2293]143  use_filtre_fft=.FALSE.
[2336]144  CALL getin('use_filtre_fft',use_filtre_fft)
[2293]145  IF(use_filtre_fft) THEN
[2336]146     WRITE(lunout,*)"FFT filter not available for sequential dynamics."
147     WRITE(lunout,*)"Your setting of variable use_filtre_fft is not used."
[2293]148  ENDIF
149
[2665]150!--- LAND MASK. THREE CASES:
[2336]151!   1) read from ocean model    file "o2a.nc"    (coupled runs)
[2665]152!   2) read from previous run   file="startphy0.nc"
153!   3) computed from topography file "Relief.nc" (masque(:,:)=-99999.)
154! In the first case, the mask from the ocean model is used compute the
[2336]155! weights to ensure ocean fractions are the same for atmosphere and ocean.
156!*******************************************************************************
[2665]157  IF(NF90_OPEN("o2a.nc", NF90_NOWRITE, nid_o2a)==NF90_NOERR) THEN
[2336]158    iret=NF90_CLOSE(nid_o2a)
159    WRITE(lunout,*)'BEWARE !! Ocean mask "o2a.nc" file found'
160    WRITE(lunout,*)'Coupled run.'
161    CALL flininfo("o2a.nc", iml_omask, jml_omask, llm_tmp, ttm_tmp, nid_o2a)
162    IF(iml_omask/=iim .OR.jml_omask/=jjp1) THEN
163      WRITE(lunout,*)'Mismatching dimensions for ocean mask'
164      WRITE(lunout,*)'iim  = ',iim ,' iml_omask = ',iml_omask
165      WRITE(lunout,*)'jjp1 = ',jjp1,' jml_omask = ',jml_omask
166      CALL abort_gcm(modname,'',1)
167    END IF
168    ALLOCATE(ocemask(iim,jjp1),lon_omask(iim,jjp1),dlon_omask(iim ))
169    ALLOCATE(ocetmp (iim,jjp1),lat_omask(iim,jjp1),dlat_omask(jjp1))
[2338]170    CALL flinopen("o2a.nc", .FALSE.,iml_omask,jml_omask,llm_tmp,               &
171                  lon_omask,lat_omask,lev,ttm_tmp,itaul,date,dt,fid)
[2336]172    CALL flinget(fid, "OceMask",    iim,jjp1,llm_tmp,ttm_tmp,1,1,ocetmp)
173    CALL flinclo(fid)
174    dlon_omask(1:iim ) = lon_omask(1:iim,1)
175    dlat_omask(1:jjp1) = lat_omask(1,1:jjp1)
176    ocemask = ocetmp
177    IF(dlat_omask(1)<dlat_omask(jml_omask)) THEN
[3815]178       DO j=1,jjp1
179          ocemask(:,j) = ocetmp(:,jjp1-j+1)
180       END DO
[2336]181    END IF
182    DEALLOCATE(ocetmp,lon_omask,lat_omask,dlon_omask,dlat_omask)
183    IF(prt_level>=1) THEN
184      WRITE(fmt,"(i4,'i1)')")iim ; fmt='('//ADJUSTL(fmt)
185      WRITE(lunout,*)'OCEAN MASK :'
186      WRITE(lunout,fmt) NINT(ocemask)
187    END IF
188    masque(1:iim,:)=1.-ocemask(:,:)
189    masque(iip1 ,:)=masque(1,:)
190    DEALLOCATE(ocemask)
[2665]191  ELSE IF(NF90_OPEN("startphy0.nc", NF90_NOWRITE, nid_sta)==NF90_NOERR) THEN
192    WRITE(lunout,*)'BEWARE !! File "startphy0.nc" found.'
193    WRITE(lunout,*)'Getting the land mask from a previous run.'
194    iret=NF90_INQ_DIMID(nid_sta,'points_physiques',nid_nph)
195    iret=NF90_INQUIRE_DIMENSION(nid_sta,nid_nph,len=nphys)
196    IF(nphys/=klon) THEN
197      WRITE(lunout,*)'Mismatching dimensions for land mask'
198      WRITE(lunout,*)'nphys  = ',nphys ,' klon = ',klon
199      iret=NF90_CLOSE(nid_sta)
200      CALL abort_gcm(modname,'',1)
201    END IF
202    ALLOCATE(masktmp(klon))
203    iret=NF90_INQ_VARID(nid_sta,'masque',nid_msk)
204    iret=NF90_GET_VAR(nid_sta,nid_msk,masktmp)
205    iret=NF90_CLOSE(nid_sta)
206    CALL gr_fi_dyn(1,klon,iip1,jjp1,masktmp,masque)
207    IF(prt_level>=1) THEN
208      WRITE(fmt,"(i4,'i1)')")iip1 ; fmt='('//ADJUSTL(fmt)
209      WRITE(lunout,*)'LAND MASK :'
210      WRITE(lunout,fmt) NINT(masque)
211    END IF
212    DEALLOCATE(masktmp)
213  ELSE
214    WRITE(lunout,*)'BEWARE !! No ocean mask "o2a.nc" file or "startphy0.nc" file found'
215    WRITE(lunout,*)'Land mask will be built from the topography file.'
216    masque(:,:)=-99999.
[2293]217  END IF
[2336]218  phis(:,:)=-99999.
[2293]219
220  IF(ok_etat0) THEN
[2336]221    WRITE(lunout,'(//)')
222    WRITE(lunout,*) '  ************************  '
223    WRITE(lunout,*) '  ***  etat0phy_netcdf ***  '
224    WRITE(lunout,*) '  ************************  '
225    CALL etat0phys_netcdf(masque,phis)
226    WRITE(lunout,'(//)')
227    WRITE(lunout,*) '  ************************  '
228    WRITE(lunout,*) '  ***  etat0dyn_netcdf ***  '
229    WRITE(lunout,*) '  ************************  '
230    CALL etat0dyn_netcdf(masque,phis)
[2293]231  END IF
[524]232
[2336]233  IF(ok_limit) THEN
234    WRITE(lunout,'(//)')
235    WRITE(lunout,*) '  *********************  '
236    WRITE(lunout,*) '  ***  Limit_netcdf ***  '
237    WRITE(lunout,*) '  *********************  '
238    WRITE(lunout,'(//)')
239    CALL limit_netcdf(masque,phis,extrap)
[2293]240  END IF
[524]241
[2336]242  WRITE(lunout,'(//)')
243  WRITE(lunout,*) '  ***************************  '
244  WRITE(lunout,*) '  ***  grilles_gcm_netcdf ***  '
245  WRITE(lunout,*) '  ***************************  '
246  WRITE(lunout,'(//)')
247  CALL grilles_gcm_netcdf_sub(masque,phis)
248
249#ifdef CPP_PARA
[1319]250  END IF
[4619]251  IF (using_xios) CALL xios_finalize
[3815]252  IF (using_mpi) call MPI_FINALIZE(ierr)
253#endif
[1319]254
255END PROGRAM ce0l
[2336]256!
257!-------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.