source: trunk/LMDZ.PLUTO/libf/phypluto/tabfi_mod.F90 @ 3811

Last change on this file since 3811 was 3772, checked in by afalco, 7 months ago

Pluto: allow to write controle in XIOS output.
AF

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