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

Last change on this file since 3300 was 2635, checked in by aslmd, 3 years ago

Introduction of cpp_mugaz_mode to decide what source should be used for cpp and mugaz. Force_cpp and check_cpp_match are deprecated (see issue la-communaut-des-mod-les-atmosph-riques-plan-taires/git-trunk#27).

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