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

Last change on this file since 2551 was 2545, checked in by romain.vande, 4 years ago

MARS Dynamico:

tab_cntrl is written by Xios output. It is now a module variable of phyetat0_mod
RV

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