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

Last change on this file since 3562 was 3552, checked in by emillour, 9 days ago

Generic PCM

Modify iostart routines to take netcdf index as an argument in order to
be used to read/write in other files than startfi.nc (preparing 1D
restart)

MM

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