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

Last change on this file since 3728 was 3728, checked in by emillour, 2 months ago

Generic PCM:
Minor fix in tabfi: the physics time step stored in the startfi file should
not overwrite the real value of dtphys (computed in and passed by the dynamics).
EM

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