source: trunk/LMDZ.GENERIC/libf/phystd/tabfi_mod.F90 @ 3922

Last change on this file since 3922 was 3908, checked in by emillour, 5 months ago

Generic PCM:
Fix some missing initializations in newstart and add
option "Lmodif=2" when calling tabi from newstart to skip some checks
which only make sense when called during a regular GCM run.
EM

File size: 23.6 KB
Line 
1MODULE tabfi_mod
2
3IMPLICIT NONE
4
5real, save :: tab_cntrl_mod(100)
6!$OMP THREADPRIVATE(tab_cntrl_mod)
7
8
9CONTAINS
10
11!=======================================================================
12      SUBROUTINE tabfi(ngrid,nid,Lmodif,tab0,day_ini,lmax,p_rad, &
13                       p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
14!=======================================================================
15!
16!   C. Hourdin 15/11/96
17!
18!   Object:        Lecture du tab_cntrl physique dans un fichier
19!   ------            et initialisation des constantes physiques
20!
21!   Arguments:
22!   ----------
23!
24!     Inputs:
25!     ------
26!
27!      - nid:    unitne logique du fichier ou on va lire le tab_cntrl   
28!                      (ouvert dans le programme appellant)
29!
30!                 si nid=0:
31!                       pas de lecture du tab_cntrl mais
32!                       Valeurs par default des constantes physiques
33!       
34!      - tab0:    Offset de tab_cntrl a partir duquel sont ranges
35!                  les parametres physiques (50 pour start_archive)
36!
37!      - Lmodif:  si on souhaite modifier les constantes  Lmodif = 1 = TRUE
38!
39!
40!     Outputs:
41!     --------
42!
43!      - day_ini: tab_cntrl(tab0+3) (Dans les cas ou l'on souhaite
44!                              comparer avec le day_ini dynamique)
45!
46!      - lmax:    tab_cntrl(tab0+2) (pour test avec nlayer)
47!
48!      - p_rad
49!      - p_omeg   !
50!      - p_g      ! Constantes physiques ayant des
51!      - p_mugaz  ! homonymes dynamiques
52!      - p_daysec !
53!
54!=======================================================================
55! to use  'getin_p'
56      use ioipsl_getin_p_mod, only: getin_p
57
58      use surfdat_h, only: emisice, iceradius, dtemisice, &
59                           emissiv
60      use comsoil_h, only: volcapa
61      use iostart, only: get_var
62      use mod_phys_lmdz_para, only: is_parallel
63      use planete_mod, only: year_day, periastr, apoastr, peri_day, &
64                             obliquit, z0
65      use comcstfi_mod, only: rad, omeg, g, mugaz, rcp, cpp, r
66      use time_phylmdz_mod, only: dtphys, daysec
67      use callkeys_mod, only: cpp_mugaz_mode
68      implicit none
69 
70      include "netcdf.inc"
71
72!-----------------------------------------------------------------------
73!   Declarations
74!-----------------------------------------------------------------------
75
76! Arguments
77! ---------
78      INTEGER,INTENT(IN) :: ngrid
79      INTEGER,INTENT(IN) :: nid ! 0: no file to read controle array from
80                                ! /= 0 : read controle array from file
81      INTEGER,INTENT(IN) :: tab0
82      INTEGER*4,INTENT(OUT) :: day_ini
83      INTEGER,INTENT(IN) :: Lmodif ! 0: do not modifications
84                                   ! 1: enable modifications
85                                   ! 2: enable modifications but no sanity test
86                                   !    (needed when called from newstart)
87      INTEGER,INTENT(OUT) :: lmax
88      REAL,INTENT(OUT) :: p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time
89
90! Variables
91! ---------
92      INTEGER,PARAMETER :: length=100
93      REAL tab_cntrl(length) ! array in which are stored the run's parameters
94      INTEGER  ierr,nvarid
95      INTEGER size
96      CHARACTER modif*20
97      LOGICAL :: found
98      CHARACTER(len=5) :: modname="tabfi"
99      real :: cntrl_daysec, cntrl_dtphys
100     
101      write(*,*)"tabfi: nid=",nid," tab0=",tab0," Lmodif=",Lmodif
102
103      IF (nid.eq.0) then
104!-----------------------------------------------------------------------
105!  Initialization of various physical constants to defaut values (nid = 0 case)
106!-----------------------------------------------------------------------
107        tab_cntrl(:)=0
108        lmax=0 ! not used anyways
109        !day_ini already set via inifis
110        time=0
111! Informations about planet for dynamics and physics
112        ! rad,cpp,g,r,rcp already initialized by inifis
113        omeg=-999.
114        call getin_p("omega",omeg)
115        if (omeg.eq.-999.) then
116          call abort_physic(modname,"Missing value for omega in def files!",1)
117        endif
118        mugaz=(8.314463/r)*1.E3
119        ! daysec already set by inifis
120        ! dtphys alread set by inifis
121! Informations about planet for the physics only
122        year_day=-999. ! length of year, in standard days
123        call getin_p("year_day",year_day)
124        if (year_day.eq.-999.) then
125          call abort_physic(modname, &
126               "Missing value for year_day in def files!",1)
127        endif
128        periastr=-999.
129        call getin_p("periastron",periastr)
130        if (periastr.eq.-999.) then
131          call abort_physic(modname, &
132               "Missing value for periastron in def files!",1)
133        endif
134        apoastr=-999.
135        call getin_p("apoastron",apoastr)
136        if (apoastr.eq.-999.) then
137          call abort_physic(modname, &
138               "Missing value for apoastron in def files!",1)
139        endif
140        peri_day=-999.
141        call getin_p("periastron_day",peri_day)
142        if (peri_day.eq.-999.) then
143          call abort_physic(modname, &
144               "Missing value for periastron date in def files!",1)
145        endif
146        obliquit=-999.
147        call getin_p("obliquity",obliquit)
148        if (obliquit.eq.-999.) then
149          call abort_physic(modname, &
150               "Missing value for obliquity in def files!",1)
151        endif
152! boundary layer and turbulence
153        z0=1.e-2 ! surface roughness length (m)
154
155! optical properties of polar caps and ground emissivity
156        emisice(:)=0
157        emissiv=0
158        iceradius(:)=1.e-6 ! mean scat radius of CO2 snow
159        dtemisice(:)=0 !time scale for snow metamorphism
160        volcapa=1000000 ! volumetric heat capacity of subsurface
161
162      ELSE
163!-----------------------------------------------------------------------
164!  Initialization of physical constants by reading array tab_cntrl(:)
165!               which contains these parameters (nid != 0 case)
166!-----------------------------------------------------------------------
167! Read 'controle' array
168!
169
170       call get_var(nid,"controle",tab_cntrl,found)
171       if (.not.found) then
172         call abort_physic(modname,"Failed reading <controle> array",1)
173       else
174         write(*,*)'tabfi: tab_cntrl',tab_cntrl
175       endif
176!
177!  Initialization of some physical constants
178! informations on physics grid
179!      if(ngrid.ne.tab_cntrl(tab0+1)) then
180!         print*,'tabfi: WARNING !!! tab_cntrl(tab0+1).ne.ngrid'
181!         print*,tab_cntrl(tab0+1),ngrid
182!      endif
183      lmax = nint(tab_cntrl(tab0+2))
184      day_ini = tab_cntrl(tab0+3)
185      time = tab_cntrl(tab0+4)
186      write (*,*) 'IN tabfi day_ini=',day_ini
187! Informations about planet for dynamics and physics
188      rad = tab_cntrl(tab0+5)
189      omeg = tab_cntrl(tab0+6)
190      g = tab_cntrl(tab0+7)
191      mugaz = tab_cntrl(tab0+8)
192      rcp = tab_cntrl(tab0+9)
193      cpp=(8.314463/(mugaz/1000.0))/rcp
194! daysec and dtphys are already set by inifis
195      cntrl_daysec = tab_cntrl(tab0+10)
196      cntrl_dtphys = tab_cntrl(tab0+11)
197! Informations about planet for the physics only
198      year_day = tab_cntrl(tab0+14)
199      periastr = tab_cntrl(tab0+15)
200      apoastr = tab_cntrl(tab0+16)
201      peri_day = tab_cntrl(tab0+17)
202      obliquit = tab_cntrl(tab0+18)
203! boundary layer and turbulence
204      z0 = tab_cntrl(tab0+19)
205! optical properties of polar caps and ground emissivity
206      emisice(1) = tab_cntrl(tab0+24)
207      emisice(2) = tab_cntrl(tab0+25)
208      emissiv    = tab_cntrl(tab0+26)
209      iceradius(1)= tab_cntrl(tab0+31) ! mean scat radius of CO2 snow (north)
210      iceradius(2)= tab_cntrl(tab0+32) ! mean scat radius of CO2 snow (south)
211      dtemisice(1)= tab_cntrl(tab0+33) !time scale for snow metamorphism (north)
212      dtemisice(2)= tab_cntrl(tab0+34) !time scale for snow metamorphism (south)
213! soil properties
214      volcapa = tab_cntrl(tab0+35) ! volumetric heat capacity
215! warn user if there has been a change in daysec or dtphys
216      if (Lmodif < 2) then
217       ! this is meaningful if called in PCM, but not when called from newstart
218       if (daysec.ne.cntrl_daysec) then
219        write(*,*) modname//" Warning: lenght of day daysec changed from ", &
220                   cntrl_daysec," to ",daysec
221       endif
222       if (dtphys.ne.cntrl_dtphys) then
223        write(*,*) modname//" Warning: time step dtphys changed from ", &
224                   cntrl_dtphys," to ",dtphys
225       endif
226      else
227        ! set dtphys and daysec to those from controle
228        daysec=cntrl_daysec
229        dtphys=cntrl_dtphys
230      endif ! of if (Lmodif < 2)
231
232!-----------------------------------------------------------------------
233!       Save some constants for later use (as routine arguments)
234!-----------------------------------------------------------------------
235      p_omeg = omeg
236      p_g = g
237      p_cpp = cpp
238      p_mugaz = mugaz
239      p_daysec = daysec
240      p_rad=rad
241
242      ENDIF    ! end of (nid = 0)
243
244!-----------------------------------------------------------------------
245!       Write physical constants to output before modifying them
246!-----------------------------------------------------------------------
247 
248   6  FORMAT(a20,e15.6,e15.6)
249   5  FORMAT(a20,f12.2,f12.2)
250 
251      write(*,*) '*****************************************************'
252      write(*,*) 'Reading tab_cntrl when calling tabfi before changes'
253      write(*,*) '*****************************************************'
254      write(*,5) '(1)        = ngrid?',tab_cntrl(tab0+1),float(ngrid)
255      write(*,5) '(2)            lmax',tab_cntrl(tab0+2),float(lmax)
256      write(*,5) '(3)         day_ini',tab_cntrl(tab0+3),float(day_ini)
257      write(*,5) '(5)             rad',tab_cntrl(tab0+5),rad
258      write(*,5) '(10)         daysec',tab_cntrl(tab0+10),daysec
259      write(*,6) '(6)            omeg',tab_cntrl(tab0+6),omeg
260      write(*,5) '(7)               g',tab_cntrl(tab0+7),g
261      write(*,5) '(8)           mugaz',tab_cntrl(tab0+8),mugaz
262      write(*,5) '(9)             rcp',tab_cntrl(tab0+9),rcp
263      write(*,6) '(11)        dtphys?',tab_cntrl(tab0+11),dtphys
264
265      write(*,5) '(14)       year_day',tab_cntrl(tab0+14),year_day
266      write(*,5) '(15)       periastr',tab_cntrl(tab0+15),periastr
267      write(*,5) '(16)        apoastr',tab_cntrl(tab0+16),apoastr
268      write(*,5) '(17)       peri_day',tab_cntrl(tab0+17),peri_day
269      write(*,5) '(18)       obliquit',tab_cntrl(tab0+18),obliquit
270
271      write(*,6) '(19)             z0',tab_cntrl(tab0+19),z0
272
273      write(*,5) '(26)        emissiv',tab_cntrl(tab0+26),emissiv
274      write(*,5) '(24)     emisice(1)',tab_cntrl(tab0+24),emisice(1)
275      write(*,5) '(25)     emisice(2)',tab_cntrl(tab0+25),emisice(2)
276      write(*,6) '(31)   iceradius(1)',tab_cntrl(tab0+31),iceradius(1)
277      write(*,6) '(32)   iceradius(2)',tab_cntrl(tab0+32),iceradius(2)
278      write(*,5) '(33)   dtemisice(1)',tab_cntrl(tab0+33),dtemisice(1)
279      write(*,5) '(34)   dtemisice(2)',tab_cntrl(tab0+34),dtemisice(2)
280
281      write(*,5) '(35)        volcapa',tab_cntrl(tab0+35),volcapa
282
283      write(*,*)
284      write(*,*) 'Lmodif in tabfi!!!!!!!',Lmodif
285
286!-----------------------------------------------------------------------
287!        Modifications...
288! NB: Modifying controls should only be done by newstart, and in seq mode
289      if ((Lmodif >= 1).and.is_parallel) then
290        write(*,*) "tabfi: Error modifying tab_control should", &
291                   " only happen in serial mode (eg: by newstart)"
292        call abort_physic("tabfi", "tab_control can be modified only by newstart or in seq mode", 1)
293      endif
294!-----------------------------------------------------------------------
295
296      IF(Lmodif >= 1) then
297
298      write(*,*)
299      write(*,*) 'Change values in tab_cntrl ? :'
300      write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
301      write(*,*) '(Current values given above)'
302      write(*,*)
303      write(*,*) '(3)          day_ini : Initial day (=0 at Ls=0)'
304      write(*,*) '(19)              z0 :  surface roughness (m)'
305      write(*,*) '(26)         emissiv : ground emissivity'
306      write(*,*) '(24 et 25)   emisice : CO2 ice max emissivity '
307      write(*,*) '(31 et 32) iceradius : mean scat radius of CO2 snow'
308      write(*,*) '(33 et 34) dtemisice : time scale for snow metamorphism'
309      write(*,*) '(35)      volcapa : soil volumetric heat capacity'
310      write(*,*) '(18)     obliquit : planet obliquity (deg)'
311      write(*,*) '(17)     peri_day : periastron date (sols since Ls=0)'
312      write(*,*) '(15)     periastr : min. star-planet dist (UA)'
313      write(*,*) '(16)     apoastr  : max. star-planet (UA)'
314      write(*,*) '(14)     year_day : length of year (in sols)'
315      write(*,*) '(5)      rad      : radius of the planet (m)'
316      write(*,*) '(6)      omeg     : planet rotation rate (rad/s)'
317      write(*,*) '(7)      g        : gravity (m/s2)'
318      write(*,*) '(8)      mugaz    : molecular mass '
319      write(*,*) '                       of the atmosphere (g/mol)'
320      write(*,*) '(9)      rcp      : r/Cp'
321      write(*,*) '(8)+(9)  calc_cpp_mugaz : r/Cp and mugaz '
322      write(*,*) '                 computed from gases.def'
323      write(*,*) '(10)     daysec   : length of a sol (s)'
324      write(*,*)
325 
326 
327      do while(modif(1:1).ne.'hello')
328        write(*,*)
329        write(*,*)
330        write(*,*) 'Changes to perform ?'
331        write(*,*) '   (enter keyword or return )'
332        write(*,*)
333        read(*,fmt='(a20)') modif
334        if (modif(1:1) .eq. ' ') goto 999
335 
336        write(*,*)
337        write(*,*) modif(1:len_trim(modif)) , ' : '
338
339        if (modif(1:len_trim(modif)) .eq. 'day_ini') then
340          write(*,*) 'current value:',day_ini
341          write(*,*) 'enter new value:'
342 101      read(*,*,iostat=ierr) day_ini
343          if(ierr.ne.0) goto 101
344          write(*,*) ' '
345          write(*,*) 'day_ini (new value):',day_ini
346
347        else if (modif(1:len_trim(modif)) .eq. 'z0') then
348          write(*,*) 'current value:',z0
349          write(*,*) 'enter new value:'
350 102      read(*,*,iostat=ierr) z0
351          if(ierr.ne.0) goto 102
352          write(*,*) ' '
353          write(*,*) ' z0 (new value):',z0
354
355        else if (modif(1:len_trim(modif)) .eq. 'emissiv') then
356          write(*,*) 'current value:',emissiv
357          write(*,*) 'enter new value:'
358 105      read(*,*,iostat=ierr) emissiv
359          if(ierr.ne.0) goto 105
360          write(*,*) ' '
361          write(*,*) ' emissiv (new value):',emissiv
362
363        else if (modif(1:len_trim(modif)) .eq. 'emisice') then
364          write(*,*) 'current value emisice(1) North:',emisice(1)
365          write(*,*) 'enter new value:'
366 106      read(*,*,iostat=ierr) emisice(1)
367          if(ierr.ne.0) goto 106
368          write(*,*)
369          write(*,*) ' emisice(1) (new value):',emisice(1)
370          write(*,*)
371
372          write(*,*) 'current value emisice(2) South:',emisice(2)
373          write(*,*) 'enter new value:'
374 107      read(*,*,iostat=ierr) emisice(2)
375          if(ierr.ne.0) goto 107
376          write(*,*)
377          write(*,*) ' emisice(2) (new value):',emisice(2)
378
379        else if (modif(1:len_trim(modif)) .eq. 'iceradius') then
380          write(*,*) 'current value iceradius(1) North:',iceradius(1)
381          write(*,*) 'enter new value:'
382 110      read(*,*,iostat=ierr) iceradius(1)
383          if(ierr.ne.0) goto 110
384          write(*,*)
385          write(*,*) ' iceradius(1) (new value):',iceradius(1)
386          write(*,*)
387
388          write(*,*) 'current value iceradius(2) South:',iceradius(2)
389          write(*,*) 'enter new value:'
390 111      read(*,*,iostat=ierr) iceradius(2)
391          if(ierr.ne.0) goto 111
392          write(*,*)
393          write(*,*) ' iceradius(2) (new value):',iceradius(2)
394
395        else if (modif(1:len_trim(modif)) .eq. 'dtemisice') then
396          write(*,*) 'current value dtemisice(1) North:',dtemisice(1)
397          write(*,*) 'enter new value:'
398 112      read(*,*,iostat=ierr) dtemisice(1)
399          if(ierr.ne.0) goto 112
400          write(*,*)
401          write(*,*) ' dtemisice(1) (new value):',dtemisice(1)
402          write(*,*)
403
404          write(*,*) 'current value dtemisice(2) South:',dtemisice(2)
405          write(*,*) 'enter new value:'
406 113      read(*,*,iostat=ierr) dtemisice(2)
407          if(ierr.ne.0) goto 113
408          write(*,*)
409          write(*,*) ' dtemisice(2) (new value):',dtemisice(2)
410
411        else if (modif(1:len_trim(modif)) .eq. 'obliquit') then
412          write(*,*) 'current value:',obliquit
413          write(*,*) 'obliquit should be 25.19 on current Mars'
414          write(*,*) 'enter new value:'
415 115      read(*,*,iostat=ierr) obliquit
416          if(ierr.ne.0) goto 115
417          write(*,*)
418          write(*,*) ' obliquit (new value):',obliquit
419
420        else if (modif(1:len_trim(modif)) .eq. 'peri_day') then
421          write(*,*) 'current value:',peri_day
422          write(*,*) 'peri_day should be 485 on current Mars'
423          write(*,*) 'enter new value:'
424 116      read(*,*,iostat=ierr) peri_day
425          if(ierr.ne.0) goto 116
426          write(*,*)
427          write(*,*) ' peri_day (new value):',peri_day
428
429        else if (modif(1:len_trim(modif)) .eq. 'periastr') then
430          write(*,*) 'current value:',periastr
431          write(*,*) 'periastr should be 1.3814 AU on present-day Mars'
432          write(*,*) 'enter new value:'
433 117      read(*,*,iostat=ierr) periastr
434          if(ierr.ne.0) goto 117
435          write(*,*)
436          write(*,*) ' periastr (new value):',periastr
437 
438        else if (modif(1:len_trim(modif)) .eq. 'apoastr') then
439          write(*,*) 'current value:',apoastr
440          write(*,*) 'apoastr should be 1.666 AU on present-day Mars'
441          write(*,*) 'enter new value:'
442 118      read(*,*,iostat=ierr) apoastr
443          if(ierr.ne.0) goto 118
444          write(*,*)
445          write(*,*) ' apoastr (new value):',apoastr
446 
447        else if (modif(1:len_trim(modif)) .eq. 'volcapa') then
448          write(*,*) 'current value:',volcapa
449          write(*,*) 'enter new value:'
450 119      read(*,*,iostat=ierr) volcapa
451          if(ierr.ne.0) goto 119
452          write(*,*)
453          write(*,*) ' volcapa (new value):',volcapa
454       
455        else if (modif(1:len_trim(modif)).eq.'rad') then
456          write(*,*) 'current value:',rad
457          write(*,*) 'enter new value:'
458 120      read(*,*,iostat=ierr) rad
459          if(ierr.ne.0) goto 120
460          write(*,*)
461          write(*,*) ' rad (new value):',rad
462
463        else if (modif(1:len_trim(modif)).eq.'omeg') then
464          write(*,*) 'current value:',omeg
465          write(*,*) 'enter new value:'
466 121      read(*,*,iostat=ierr) omeg
467          if(ierr.ne.0) goto 121
468          write(*,*)
469          write(*,*) ' omeg (new value):',omeg
470       
471        else if (modif(1:len_trim(modif)).eq.'g') then
472          write(*,*) 'current value:',g
473          write(*,*) 'enter new value:'
474 122      read(*,*,iostat=ierr) g
475          if(ierr.ne.0) goto 122
476          write(*,*)
477          write(*,*) ' g (new value):',g
478
479        else if (modif(1:len_trim(modif)).eq.'mugaz') then
480          write(*,*) 'current value:',mugaz
481          write(*,*) 'enter new value:'
482 123      read(*,*,iostat=ierr) mugaz
483          if(ierr.ne.0) goto 123
484          write(*,*)
485          write(*,*) ' mugaz (new value):',mugaz
486          r=8.314463/(mugaz/1000.0)
487          write(*,*) ' R (new value):',r
488
489        else if (modif(1:len_trim(modif)).eq.'rcp') then
490          write(*,*) 'current value:',rcp
491          write(*,*) 'enter new value:'
492 124      read(*,*,iostat=ierr) rcp
493          if(ierr.ne.0) goto 124
494          write(*,*)
495          write(*,*) ' rcp (new value):',rcp
496          r=8.314463/(mugaz/1000.0)
497          cpp=r/rcp
498          write(*,*) ' cpp (new value):',cpp
499
500        else if (modif(1:len_trim(modif)).eq.'calc_cpp_mugaz') then
501          write(*,*) 'current value rcp, mugaz:',rcp,mugaz
502          cpp_mugaz_mode = 2
503          call su_gases
504          call calc_cpp_mugaz
505          write(*,*)
506          write(*,*) ' cpp (new value):',cpp
507          write(*,*) ' mugaz (new value):',mugaz
508          r=8.314463/(mugaz/1000.0)
509          rcp=r/cpp
510          write(*,*) ' rcp (new value):',rcp
511         
512        else if (modif(1:len_trim(modif)).eq.'daysec') then
513          write(*,*) 'current value:',daysec
514          write(*,*) 'enter new value:'
515 125      read(*,*,iostat=ierr) daysec
516          if(ierr.ne.0) goto 125
517          write(*,*)
518          write(*,*) ' daysec (new value):',daysec
519
520!         added by RW!
521        else if (modif(1:len_trim(modif)).eq.'year_day') then
522          write(*,*) 'current value:',year_day
523          write(*,*) 'enter new value:'
524 126      read(*,*,iostat=ierr) year_day
525          if(ierr.ne.0) goto 126
526          write(*,*)
527          write(*,*) ' year_day (new value):',year_day
528
529        endif
530      enddo ! of do while(modif(1:1).ne.'hello')
531
532 999  continue
533
534!----------------------------------------------------------------------
535!       Write values of physical constants after modifications
536!-----------------------------------------------------------------------
537 
538      write(*,*) '*****************************************************'
539      write(*,*) 'Reading tab_cntrl when calling tabfi AFTER changes'
540      write(*,*) '*****************************************************'
541      write(*,5) '(1)        = ngrid?',tab_cntrl(tab0+1),float(ngrid)
542      write(*,5) '(2)            lmax',tab_cntrl(tab0+2),float(lmax)
543      write(*,5) '(3)         day_ini',tab_cntrl(tab0+3),float(day_ini)
544      write(*,5) '(5)             rad',tab_cntrl(tab0+5),rad
545      write(*,5) '(10)         daysec',tab_cntrl(tab0+10),daysec
546      write(*,6) '(6)            omeg',tab_cntrl(tab0+6),omeg
547      write(*,5) '(7)               g',tab_cntrl(tab0+7),g
548      write(*,5) '(8)           mugaz',tab_cntrl(tab0+8),mugaz
549      write(*,5) '(9)             rcp',tab_cntrl(tab0+9),rcp
550      write(*,6) '(11)        dtphys?',tab_cntrl(tab0+11),dtphys
551 
552      write(*,5) '(14)       year_day',tab_cntrl(tab0+14),year_day
553      write(*,5) '(15)       periastr',tab_cntrl(tab0+15),periastr
554      write(*,5) '(16)        apoastr',tab_cntrl(tab0+16),apoastr
555      write(*,5) '(17)       peri_day',tab_cntrl(tab0+17),peri_day
556      write(*,5) '(18)       obliquit',tab_cntrl(tab0+18),obliquit
557 
558      write(*,6) '(19)             z0',tab_cntrl(tab0+19),z0
559 
560      write(*,5) '(26)        emissiv',tab_cntrl(tab0+26),emissiv
561      write(*,5) '(24)     emisice(1)',tab_cntrl(tab0+24),emisice(1)
562      write(*,5) '(25)     emisice(2)',tab_cntrl(tab0+25),emisice(2)
563      write(*,6) '(31)   iceradius(1)',tab_cntrl(tab0+31),iceradius(1)
564      write(*,6) '(32)   iceradius(2)',tab_cntrl(tab0+32),iceradius(2)
565      write(*,5) '(33)   dtemisice(1)',tab_cntrl(tab0+33),dtemisice(1)
566      write(*,5) '(34)   dtemisice(2)',tab_cntrl(tab0+34),dtemisice(2)
567 
568      write(*,5) '(35)        volcapa',tab_cntrl(tab0+35),volcapa
569
570      write(*,*) 
571      write(*,*)
572
573      ENDIF ! of if (Lmodif >= 1)
574
575!-----------------------------------------------------------------------
576!       Save some constants for later use (as routine arguments)
577!-----------------------------------------------------------------------
578      p_omeg = omeg
579      p_g = g
580      p_cpp = cpp
581      p_mugaz = mugaz
582      p_daysec = daysec
583      p_rad=rad
584
585
586      END SUBROUTINE tabfi
587
588      subroutine ini_tab_controle_dyn_xios(idayref)
589      ! Initialize controle variable for XIOS & diagfi
590
591        use comcstfi_mod,        only: g, mugaz, omeg, rad, rcp
592        use time_phylmdz_mod,  only: daysec, dtphys
593        use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev
594
595        implicit none
596
597        integer*4, intent(in) :: idayref ! date (initial date for this run)
598
599        integer :: length, l
600        parameter (length = 100)
601        real    :: tab_cntrl(length) ! run parameters are stored in this array
602
603        do l = 1,length
604            tab_cntrl(l) = 0.
605        enddo
606        tab_cntrl(1)  = real(nbp_lon)
607        tab_cntrl(2)  = real(nbp_lat-1)
608        tab_cntrl(3)  = real(nbp_lev)
609        tab_cntrl(4)  = real(idayref)
610        tab_cntrl(5)  = rad
611        tab_cntrl(6)  = omeg
612        tab_cntrl(7)  = g
613        tab_cntrl(8)  = mugaz
614        tab_cntrl(9)  = rcp
615        tab_cntrl(10) = daysec
616        tab_cntrl(11) = dtphys
617
618        tab_cntrl_mod = tab_cntrl
619
620      end subroutine ini_tab_controle_dyn_xios
621
622end module tabfi_mod
Note: See TracBrowser for help on using the repository browser.