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

Last change on this file since 2616 was 2578, checked in by romain.vande, 3 years ago

First stage of implementing Open_MP in the physic.
So far it can initialyse physic and run with all routines at .FALSE.

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