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

Last change on this file since 3817 was 3773, checked in by afalco, 6 weeks ago

Generic: allow to write controle in XIOS output.
AF

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