source: trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90 @ 2566

Last change on this file since 2566 was 2562, checked in by cmathe, 4 years ago

GCM MARS: CO2 clouds microphysics improvements

File size: 20.0 KB
RevLine 
[1944]1module phyetat0_mod
2
3implicit none
[2545]4  real,save :: tab_cntrl_mod(100)
[1944]5
6contains
7
[1130]8subroutine phyetat0 (fichnom,tab0,Lmodif,nsoil,ngrid,nlay,nq, &
[1944]9                     day_ini,time0,tsurf,tsoil,albedo,emis,q2,qsurf,co2ice, &
[2562]10                     tauscaling,totcloudfrac,wstar,watercap)
[1922]11
[1621]12  use tracer_mod, only: noms ! tracer names
[1130]13  use surfdat_h, only: phisfi, albedodat, z0, z0_default,&
[2079]14                       zmea, zstd, zsig, zgam, zthe, hmons, summit, base
[1130]15  use iostart, only: nid_start, open_startphy, close_startphy, &
16                     get_field, get_var, inquire_field, &
17                     inquire_dimension, inquire_dimension_length
[2220]18  use nonoro_gwd_ran_mod, only: du_nonoro_gwd, dv_nonoro_gwd
[2265]19  use compute_dtau_mod, only: dtau
[2417]20  use dust_rad_adjust_mod, only: dust_rad_adjust_prev,dust_rad_adjust_next
21  use dust_param_mod, only: dustscaling_mode
[2281]22  USE ioipsl_getin_p_mod, ONLY : getin_p
23
[1130]24  implicit none
[1974]25 
26  include "callkeys.h"
[1130]27!======================================================================
28! Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
[2562]29!  Adaptation � Mars : Yann Wanherdrick
[1130]30! Objet: Lecture de l etat initial pour la physique
31! Modifs: Aug.2010 EM : use NetCDF90 to load variables (enables using
32!                      r4 or r8 restarts independently of having compiled
33!                      the GCM in r4 or r8)
34!         June 2013 TN : Possibility to read files with a time axis
35!         November 2013 EM : Enabeling parallel, using iostart module
[2281]36!         March 2020 AD: Enabling initialization of physics without startfi
37!                        flag: startphy_file
[1130]38!======================================================================
39  INTEGER nbsrf !Mars nbsrf a 1 au lieu de 4
40  PARAMETER (nbsrf=1) ! nombre de sous-fractions pour une maille
41!======================================================================
42!  Arguments:
43!  ---------
44!  inputs:
[2281]45!  logical,intent(in) :: startphy_file ! .true. if reading start file
[1130]46  character*(*),intent(in) :: fichnom ! "startfi.nc" file
47  integer,intent(in) :: tab0
48  integer,intent(in) :: Lmodif
49  integer,intent(in) :: nsoil ! # of soil layers
50  integer,intent(in) :: ngrid ! # of atmospheric columns
51  integer,intent(in) :: nlay ! # of atmospheric layers
52  integer,intent(in) :: nq
53  integer :: day_ini
54  real :: time0
55
56!  outputs:
57  real,intent(out) :: tsurf(ngrid) ! surface temperature
58  real,intent(out) :: tsoil(ngrid,nsoil) ! soil temperature
[1944]59  real,intent(out) :: albedo(ngrid,2) ! surface albedo
[1130]60  real,intent(out) :: emis(ngrid) ! surface emissivity
61  real,intent(out) :: q2(ngrid,nlay+1) !
62  real,intent(out) :: qsurf(ngrid,nq) ! tracers on surface
63  real,intent(out) :: co2ice(ngrid) ! co2 ice cover
[1208]64  real,intent(out) :: tauscaling(ngrid) ! dust conversion factor
[1711]65  real,intent(out) :: totcloudfrac(ngrid) ! total cloud fraction
[1944]66  real,intent(out) :: wstar(ngrid) ! Max vertical velocity in thermals (m/s)
[2260]67  real,intent(out) :: watercap(ngrid) ! h2o_ice_cover
[1130]68!======================================================================
69!  Local variables:
70
71      real surffield(ngrid) ! to temporarily store a surface field
72      real xmin,xmax ! to display min and max of a field
73!
74      INTEGER ig,iq,lmax
75      INTEGER nid, nvarid
76      INTEGER ierr, i, nsrf
77!      integer isoil
78!      INTEGER length
79!      PARAMETER (length=100)
80      CHARACTER*7 str7
81      CHARACTER*2 str2
82      CHARACTER*1 yes
83!
84      REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec
85      INTEGER nqold
86
87! flag which identifies if 'startfi.nc' file is using old names (qsurf01,...)
88      logical :: oldtracernames=.false.
89      integer :: count
90      character(len=30) :: txt ! to store some text
91
92! specific for time
93      REAL,ALLOCATABLE :: time(:) ! times stored in start
94      INTEGER timelen ! number of times stored in the file
95      INTEGER indextime ! index of selected time
96
97      INTEGER :: edges(3),corner(3)
98      LOGICAL :: found
99
[1525]100      REAL :: timestart ! to pick which initial state to start from
[2281]101      REAL :: surfemis  ! constant emissivity when no startfi
102      REAL :: surfalbedo  ! constant emissivity when no startfi
103      CHARACTER(len=5) :: modname="phyetat0"
[1525]104
[2281]105write(*,*) "phyetat0: startphy_file", startphy_file
[1130]106
[2281]107if (startphy_file) then
108   ! open physics initial state file:
109   call open_startphy(fichnom)
110   ! possibility to modify tab_cntrl in tabfi
111   write(*,*)
112   write(*,*) 'TABFI in phyeta0: Lmodif=',Lmodif," tab0=",tab0
113   call tabfi (nid_start,Lmodif,tab0,day_ini,lmax,p_rad, &
114               p_omeg,p_g,p_mugaz,p_daysec,time0)
115else ! "academic" initialization of planetary parameters via tabfi
116   call tabfi (0,0,0,day_ini,lmax,p_rad, &
117               p_omeg,p_g,p_mugaz,p_daysec,time0)
118endif ! of if (startphy_file)
[1130]119
[2281]120if (startphy_file) then
121   ! Load surface geopotential:
122   call get_field("phisfi",phisfi,found)
123   if (.not.found) then
124     call abort_physic(modname, &
125                "phyetat0: Failed loading <phisfi>",1)
126   endif
[1130]127else
[2281]128  phisfi(:)=0.
129endif ! of if (startphy_file)
130write(*,*) "phyetat0: surface geopotential <phisfi> range:", &
131               minval(phisfi), maxval(phisfi)
[1130]132
133
[2281]134if (startphy_file) then
135   ! Load bare ground albedo:
136   call get_field("albedodat",albedodat,found)
137   if (.not.found) then
138     call abort_physic(modname, &
139                "phyetat0: Failed loading <albedodat>",1)
140   endif
141else ! If no startfi file, use parameter surfalbedo in def file
[2294]142  surfalbedo=0.1
143  call getin_p("surfalbedo_without_startfi",surfalbedo)
144  print*,"surfalbedo_without_startfi",surfalbedo
[2281]145  albedodat(:)=surfalbedo
146endif ! of if (startphy_file)
147write(*,*) "phyetat0: Bare ground albedo <albedodat> range:", &
[1130]148             minval(albedodat), maxval(albedodat)
149
150! ZMEA
[2281]151if (startphy_file) then
152   call get_field("ZMEA",zmea,found)
153   if (.not.found) then
154     call abort_physic(modname, &
155                "phyetat0: Failed loading <ZMEA>",1)
156   endif
[1130]157else
[2281]158  zmea(:)=0.
159endif ! of if (startphy_file)
160write(*,*) "phyetat0: <ZMEA> range:", &
161            minval(zmea), maxval(zmea)
[1130]162
163
164! ZSTD
[2281]165if (startphy_file) then
166   call get_field("ZSTD",zstd,found)
167   if (.not.found) then
168     call abort_physic(modname, &
169                "phyetat0: Failed loading <ZSTD>",1)
170   endif
[1130]171else
[2281]172  zstd(:)=0.
173endif ! of if (startphy_file)
174write(*,*) "phyetat0: <ZSTD> range:", &
175            minval(zstd), maxval(zstd)
[1130]176
177
178! ZSIG
[2281]179if (startphy_file) then
180   call get_field("ZSIG",zsig,found)
181   if (.not.found) then
182     call abort_physic(modname, &
183                "phyetat0: Failed loading <ZSIG>",1)
184   endif
[1130]185else
[2281]186  zsig(:)=0.
187endif ! of if (startphy_file)
188write(*,*) "phyetat0: <ZSIG> range:", &
189            minval(zsig), maxval(zsig)
[1130]190
191
192! ZGAM
[2281]193if (startphy_file) then
194   call get_field("ZGAM",zgam,found)
195   if (.not.found) then
196     call abort_physic(modname, &
197                "phyetat0: Failed loading <ZGAM>",1)
198   endif
[1130]199else
[2281]200  zgam(:)=0.
201endif ! of if (startphy_file)
202write(*,*) "phyetat0: <ZGAM> range:", &
203            minval(zgam), maxval(zgam)
[1130]204
205
206! ZTHE
[2281]207if (startphy_file) then
208   call get_field("ZTHE",zthe,found)
209   if (.not.found) then
210     call abort_physic(modname, &
211                "phyetat0: Failed loading <ZTHE>",1)
212   endif
[1130]213else
[2281]214  zthe(:)=0.
215endif ! of if (startphy_file)
216write(*,*) "phyetat0: <ZTHE> range:", &
[1130]217             minval(zthe), maxval(zthe)
218
[2281]219! hmons
220if (startphy_file) then
221   call get_field("hmons",hmons,found)
222   if (.not.found) then
223     write(*,*) "WARNING: phyetat0: Failed loading <hmons>"
224     if (rdstorm) then
225     call abort_physic(modname, &
226                "phyetat0: Failed loading <hmons>",1)
227     else
228       write(*,*) "will continue anyway..."
229       write(*,*) "because you may not need it."
230       hmons(:)=0.
231     end if ! if (rdstorm)
232   else
233     do ig=1,ngrid
234       if (hmons(ig) .eq. -999999.)  hmons(ig)=0.
235     enddo
236   endif ! (.not.found)
[1974]237else
[2281]238   hmons(:)=0.
239endif ! if (startphy_file)
240write(*,*) "phyetat0: <hmons> range:", &
241            minval(hmons), maxval(hmons)
[2079]242
[2281]243
244! summit
245if (startphy_file) then
246   call get_field("summit",summit,found)
247   if (.not.found) then
248     write(*,*) "WARNING: phyetat0: Failed loading <summit>"
249     if (rdstorm) then
250     call abort_physic(modname, &
251                "phyetat0: Failed loading <summit>",1)
252     else
253       write(*,*) "will continue anyway..."
254       write(*,*) "because you may not need it."
255       summit(:)=0.
256     end if
257   else
258     do ig=1,ngrid
259       if (summit(ig) .eq. -999999.)  summit(ig)=0.
260     enddo
261   endif ! if (.not.found)
[2079]262else
[2281]263   summit(:)=0. 
264endif ! if (startphy_file)
265write(*,*) "phyetat0: <summit> range:", &
266            minval(summit), maxval(summit)
[2079]267
[2281]268
269! base
270if (startphy_file) then
271   call get_field("base",base,found)
272   if (.not.found) then
273     write(*,*) "WARNING: phyetat0: Failed loading <base>"
274     if (rdstorm) then
275     call abort_physic(modname, &
276                "phyetat0: Failed loading <base>",1)
277     else
278       write(*,*) "will continue anyway..."
279       write(*,*) "because you may not need it."
280       base(:)=0.
281     end if
282   else
283     do ig=1,ngrid
284       if (base(ig) .eq. -999999.)  base(ig)=0.
285     enddo
286   endif ! if(.not.found)
[2079]287else
[2281]288   base(:)=0.
289endif ! if (startphy_file)
290write(*,*) "phyetat0: <base> range:", &
291            minval(base), maxval(base)
292
[1130]293! Time axis
[1525]294! obtain timestart from run.def
295timestart=-9999 ! default value
[2281]296call getin_p("timestart",timestart)
297if (startphy_file) then
298   found=inquire_dimension("Time")
299   if (.not.found) then
300     indextime = 1
301     write(*,*) "phyetat0: No time axis found in "//trim(fichnom)
302   else
303     write(*,*) "phyetat0: Time axis found in "//trim(fichnom)
304     timelen=inquire_dimension_length("Time")
305     allocate(time(timelen))
306     ! load "Time" array:
307     call get_var("Time",time,found)
308     if (.not.found) then
309     call abort_physic(modname, &
310                "phyetat0: Failed loading <Time>",1)
311     endif
312     ! seclect the desired time index
313     IF (timestart .lt. 0) THEN  ! default: we use the last time value
314       indextime = timelen
315     ELSE  ! else we look for the desired value in the time axis
316       indextime = 0
317       DO i=1,timelen
318         IF (abs(time(i) - timestart) .lt. 0.01) THEN
319           indextime = i
320           EXIT
321         ENDIF
322       ENDDO
323       IF (indextime .eq. 0) THEN
324         PRINT*, "Time", timestart," is not in "//trim(fichnom)//"!!"
325         PRINT*, "Stored times are:"
326         DO i=1,timelen
327            PRINT*, time(i)
328         ENDDO
329         call abort_physic(modname,"phyetat0: Time error",1)
330       ENDIF
331     ENDIF ! of IF (timestart .lt. 0)
332     ! In startfi the absolute date is day_ini + time0 + time
333     ! For now on, in the GCM physics, it is day_ini + time0
334     time0 = time(indextime) + time0
335     day_ini = day_ini + INT(time0)
336     time0 = time0 - INT(time0)   
337     PRINT*, "phyetat0: Selected time ",time(indextime), &
338             " at index ",indextime
339     DEALLOCATE(time)
340   endif ! of if Time not found in file
[2545]341
342   call ini_tab_controle_dyn_xios(day_ini)
343
[2281]344else
[1130]345  indextime = 1
[2281]346endif ! if (startphy_file)
[1130]347
348! CO2 ice cover
[2281]349if (startphy_file) then
350   call get_field("co2ice",co2ice,found,indextime)
351   if (.not.found) then
352     call abort_physic(modname, &
353                "phyetat0: Failed loading <co2ice>",1)
354   endif
[1130]355else
[2281]356   co2ice(:)=0.
357endif !if (startphy_file)
358write(*,*) "phyetat0: CO2 ice cover <co2ice> range:", &
359            minval(co2ice), maxval(co2ice)
[1130]360
[1208]361! Dust conversion factor
[2281]362if (startphy_file) then
363   call get_field("tauscaling",tauscaling,found,indextime)
364   if (.not.found) then
365     write(*,*) "phyetat0: <tauscaling> not in file"
366     tauscaling(:) = 1
367   endif
[1208]368else
[2281]369   tauscaling(:) = 1
370endif ! if (startphy_file)
371write(*,*) "phyetat0: dust conversion factor <tauscaling> range:", &
372            minval(tauscaling), maxval(tauscaling)
[1208]373
[2417]374! dust_rad_adjust_* for radiative rescaling of dust
375if (dustscaling_mode==2) then
376 if (startphy_file) then
377   call get_field("dust_rad_adjust_prev",dust_rad_adjust_prev,found,indextime)
378   if (.not.found) then
379     write(*,*) "phyetat0: <dust_rad_adjust_prev> not in file; set to 1"
380     dust_rad_adjust_prev(:) = 1
381   endif
382   call get_field("dust_rad_adjust_next",dust_rad_adjust_next,found,indextime)
383   if (.not.found) then
384     write(*,*) "phyetat0: <dust_rad_adjust_next> not in file; set to 1"
385     dust_rad_adjust_next(:) = 1
386   endif
387 else
388   dust_rad_adjust_prev(:)= 0
389   dust_rad_adjust_next(:)= 0
390 endif ! if (startphy_file)
391 write(*,*) "phyetat0: radiative scaling coeff <dust_rad_adjust_prev> range:", &
392            minval(dust_rad_adjust_prev), maxval(dust_rad_adjust_prev)
393 write(*,*) "phyetat0: radiative scaling coeff <dust_rad_adjust_next> range:", &
394            minval(dust_rad_adjust_next), maxval(dust_rad_adjust_next)
395endif ! of if (dustscaling_mode==2)
[2281]396
[2265]397! dtau: opacity difference between GCM and dust scenario
[2281]398if (startphy_file) then
399   call get_field("dtau",dtau,found,indextime)
400   if (.not.found) then
401     write(*,*) "phyetat0: <dtau> not in file; set to zero"
402     dtau(:) = 0
403   endif
[2265]404else
[2281]405   dtau(:)= 0
406endif ! if (startphy_file)
407write(*,*) "phyetat0: opacity diff wrt scenario <dtau> range:", &
408            minval(dtau), maxval(dtau)
[2265]409
410
[1711]411! Sub-grid cloud fraction
[2281]412if (startphy_file) then
413   call get_field("totcloudfrac",totcloudfrac,found,indextime)
414   if (.not.found) then
415     write(*,*) "phyetat0: <totcloudfrac> not in file WARNING put to 1"
416     totcloudfrac(:) = 1.0 !valeur par defaut (CLFfixval par defaut)
417   endif
[1711]418else
[2281]419   totcloudfrac(:)=1.0
420endif ! if (startphy_file)
421write(*,*) "phyetat0: total cloud fraction <totcloudfrac> range:", &
422            minval(totcloudfrac), maxval(totcloudfrac)
[1208]423
[2281]424
[1944]425! Max vertical velocity in thermals
[2281]426if (startphy_file) then
427   call get_field("wstar",wstar,found,indextime)
428   if (.not.found) then
429     write(*,*) "phyetat0: <wstar> not in file! Set to zero"
430     wstar(:)=0
431   endif
[1944]432else
[2281]433   wstar(:)=0
434endif ! if (startphy_file)
435write(*,*) "phyetat0: Max vertical velocity in thermals <wstar> range:", &
436            minval(wstar),maxval(wstar)
[1944]437
[2281]438
[1130]439! Surface temperature :
[2281]440if (startphy_file) then !tsurf
441   call get_field("tsurf",tsurf,found,indextime)
442   if (.not.found) then
443     call abort_physic(modname, &
444                "phyetat0: Failed loading <tsurf>",1)
445   endif
[1130]446else
[2281]447  tsurf(:)=0. ! will be updated afterwards in physiq !
448endif ! of if (startphy_file)
449write(*,*) "phyetat0: Surface temperature <tsurf> range:", &
450            minval(tsurf), maxval(tsurf)
[1130]451
[1944]452! Surface albedo
[2281]453if (startphy_file) then
454   call get_field("albedo",albedo(:,1),found,indextime)
455   if (.not.found) then
456     write(*,*) "phyetat0: Failed loading <albedo>"
457     albedo(:,1)=albedodat(:)
458   endif
[1944]459else
[2281]460   albedo(:,1)=albedodat(:)
461endif ! of if (startphy_file)
462write(*,*) "phyetat0: Surface albedo <albedo> range:", &
463            minval(albedo(:,1)), maxval(albedo(:,1))
[1944]464albedo(:,2)=albedo(:,1)
465
[1130]466! Surface emissivity
[2281]467if (startphy_file) then
468   call get_field("emis",emis,found,indextime)
469   if (.not.found) then
470     call abort_physic(modname, &
471                "phyetat0: Failed loading <emis>",1)
472   endif
[1130]473else
[2281]474  ! If no startfi file, use parameter surfemis in def file
[2294]475  surfemis=0.95
476  call getin_p("surfemis_without_startfi",surfemis)
477  print*,"surfemis_without_startfi",surfemis
[2281]478  emis(:)=surfemis
479endif ! of if (startphy_file)
480write(*,*) "phyetat0: Surface emissivity <emis> range:", &
481            minval(emis), maxval(emis)
[1130]482
483
484! surface roughness length (NB: z0 is a common in surfdat_h)
[2281]485if (startphy_file) then
486   call get_field("z0",z0,found)
487   if (.not.found) then
488     write(*,*) "phyetat0: Failed loading <z0>"
489     write(*,*) 'will use constant value of z0_default:',z0_default
490     z0(:)=z0_default
491   endif
[1130]492else
[2281]493   z0(:)=z0_default
494endif ! of if (startphy_file)
495write(*,*) "phyetat0: Surface roughness <z0> range:", &
496            minval(z0), maxval(z0)
[1130]497
498
499! pbl wind variance
[2281]500if (startphy_file) then
501   call get_field("q2",q2,found,indextime)
502   if (.not.found) then
503     call abort_physic(modname, &
504                "phyetat0: Failed loading <q2>",1)
505   endif
[1130]506else
[2281]507  q2(:,:)=0.
508endif ! of if (startphy_file)
509write(*,*) "phyetat0: PBL wind variance <q2> range:", &
510            minval(q2), maxval(q2)
[1130]511
[2281]512! Non-orographic gravity waves
513if (startphy_file) then
514   call get_field("du_nonoro_gwd",du_nonoro_gwd,found,indextime)
515   if (.not.found) then
516      write(*,*) "phyetat0: <du_nonoro_gwd> not in file"
517      du_nonoro_gwd(:,:)=0.
518   endif
[2220]519else
[2281]520du_nonoro_gwd(:,:)=0.
521endif ! if (startphy_file)
522write(*,*) "phyetat0: Memory of zonal wind tendency due to non-orographic GW"
523write(*,*) " <du_nonoro_gwd> range:", &
524             minval(du_nonoro_gwd), maxval(du_nonoro_gwd)
[1130]525
[2281]526if (startphy_file) then
527   call get_field("dv_nonoro_gwd",dv_nonoro_gwd,found,indextime)
528   if (.not.found) then
529      write(*,*) "phyetat0: <dv_nonoro_gwd> not in file"
530      dv_nonoro_gwd(:,:)=0.
531   endif
532else ! ! if (startphy_file)
533dv_nonoro_gwd(:,:)=0.
534endif ! if (startphy_file)
535write(*,*) "phyetat0: Memory of meridional wind tendency due to non-orographic GW"
536write(*,*) " <dv_nonoro_gwd> range:", &
537             minval(dv_nonoro_gwd), maxval(dv_nonoro_gwd)
[2220]538
[1130]539! tracer on surface
540if (nq.ge.1) then
541  do iq=1,nq
[1621]542    txt=noms(iq)
[1130]543    if (txt.eq."h2o_vap") then
544      ! There is no surface tracer for h2o_vap;
545      ! "h2o_ice" should be loaded instead
546      txt="h2o_ice"
547      write(*,*) 'phyetat0: loading surface tracer', &
548                           ' h2o_ice instead of h2o_vap'
[2260]549      write(*,*) 'iq = ', iq
[1130]550    endif
[2312]551
552    if (hdo) then
553    if (txt.eq."hdo_vap") then
554      txt="hdo_ice"
555      write(*,*) 'phyetat0: loading surface tracer', &
556                           ' hdo_ice instead of hdo_vap'
557    endif
558    endif !hdo
559
[2281]560    if (startphy_file) then
561       call get_field(txt,qsurf(:,iq),found,indextime)
562       if (.not.found) then
563         write(*,*) "phyetat0: Failed loading <",trim(txt),">"
564         write(*,*) "         ",trim(txt)," is set to zero"
565       endif
[1130]566    else
[2281]567      qsurf(:,iq)=0.
568    endif ! of if (startphy_file)
569    write(*,*) "phyetat0: Surface tracer <",trim(txt),"> range:", &
[1130]570                 minval(qsurf(:,iq)), maxval(qsurf(:,iq))
[2281]571  enddo ! of do iq=1,nq
[1130]572endif ! of if (nq.ge.1)
573
[2281]574if (startphy_file) then
575   call get_field("watercap",watercap,found,indextime)
576   if (.not.found) then
577     write(*,*) "phyetat0: Failed loading <watercap> : ", &
578                          "<watercap> is set to zero"
579     watercap(:)=0
[2260]580
[2281]581     write(*,*) 'Now transfer negative surface water ice to', &
582                ' watercap !'
583     if (nq.ge.1) then
584      do iq=1,nq
585       txt=noms(iq)
586       if (txt.eq."h2o_ice") then
587         do ig=1,ngrid
588          if (qsurf(ig,iq).lt.0.0) then
589             watercap(ig) = qsurf(ig,iq)
590             qsurf(ig,iq) = 0.0
591          end if
592         end do
593       endif
[2378]594
595       if (txt.eq."hdo_ice") then
596         do ig=1,ngrid
597          if (qsurf(ig,iq).lt.0.0) then
598             qsurf(ig,iq) = 0.0
599          end if
600         end do
601       endif
602
[2281]603      enddo
604     endif ! of if (nq.ge.1)
605   endif ! of if (.not.found)
[2264]606else
[2281]607   watercap(:)=0
608endif ! of if (startphy_file)
609write(*,*) "phyetat0: Surface water ice <watercap> range:", &
610            minval(watercap), maxval(watercap)
611 
[2260]612
613
[2281]614if (startphy_file) then
615  ! Call to soil_settings, in order to read soil temperatures,
616  ! as well as thermal inertia and volumetric heat capacity
617  call soil_settings(nid_start,ngrid,nsoil,tsurf,tsoil,indextime)
618endif ! of if (startphy_file)
[1130]619!
620! close file:
621!
[2281]622if (startphy_file) call close_startphy
[1130]623
[1944]624end subroutine phyetat0
625
[2545]626
627subroutine ini_tab_controle_dyn_xios(idayref)
628
629      USE comcstfi_h, only: g, mugaz, omeg, rad, rcp
630      USE time_phylmdz_mod, ONLY: hour_ini, daysec, dtphys
631      USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, nbp_lev
632      IMPLICIT NONE
633
634
635      INTEGER*4,intent(in) :: idayref  ! date (initial date for this run)
636
637
638      INTEGER length,l
639      parameter (length = 100)
640      REAL tab_cntrl(length) ! run parameters are stored in this array
641
642      DO l=1,length
643         tab_cntrl(l)=0.
644      ENDDO
645      tab_cntrl(1)  = real(nbp_lon)
646      tab_cntrl(2)  = real(nbp_lat-1)
647      tab_cntrl(3)  = real(nbp_lev)
648      tab_cntrl(4)  = real(idayref)
649      tab_cntrl(5)  = rad
650      tab_cntrl(6)  = omeg
651      tab_cntrl(7)  = g
652      tab_cntrl(8)  = mugaz
653      tab_cntrl(9)  = rcp
654      tab_cntrl(10) = daysec
655      tab_cntrl(11) = dtphys
656
657      tab_cntrl(27) = hour_ini
658
659      tab_cntrl_mod = tab_cntrl
660
661end subroutine ini_tab_controle_dyn_xios
662
[1944]663end module phyetat0_mod
Note: See TracBrowser for help on using the repository browser.