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

Last change on this file since 3477 was 3477, checked in by afalco, 4 weeks ago

Pluto PCM: fixed some variables imports for 'libphy' compilation (to compile with DYNAMICO).
AF

File size: 22.5 KB
RevLine 
[3184]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!
[3455]14!   Object:        Lecture du tab_cntrl physique dans un fichier
[3184]15!   ------            et initialisation des constantes physiques
16!
17!   Arguments:
18!   ----------
19!
20!     Inputs:
21!     ------
22!
[3455]23!      - nid:    unitne logique du fichier ou on va lire le tab_cntrl
24!                      (ouvert dans le programme appellant)
[3184]25!
26!                 si nid=0:
27!                       pas de lecture du tab_cntrl mais
28!                       Valeurs par default des constantes physiques
[3455]29!
30!      - tab0:    Offset de tab_cntrl a partir duquel sont ranges
[3184]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   !
[3455]46!      - p_g      ! Constantes physiques ayant des
[3184]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, &
[3455]60                             obliquit, z0, lmixmin, emin_turb, &
61                             tpal, adjust
[3184]62      use comcstfi_mod, only: rad, omeg, g, mugaz, rcp, cpp, r
63      use time_phylmdz_mod, only: dtphys, daysec
64      use callkeys_mod, only: cpp_mugaz_mode
65      implicit none
[3455]66
[3184]67      include "netcdf.inc"
68
69!-----------------------------------------------------------------------
70!   Declarations
71!-----------------------------------------------------------------------
72
73! Arguments
74! ---------
75      INTEGER,INTENT(IN) :: ngrid,nid,tab0
76      INTEGER*4,INTENT(OUT) :: day_ini
77      INTEGER,INTENT(IN) :: Lmodif
78      INTEGER,INTENT(OUT) :: lmax
79      REAL,INTENT(OUT) :: p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time
80
81! Variables
82! ---------
83      INTEGER,PARAMETER :: length=100
84      REAL tab_cntrl(length) ! array in which are stored the run's parameters
85      INTEGER  ierr,nvarid
86      INTEGER size
87      CHARACTER modif*20
88      LOGICAL :: found
89      CHARACTER(len=5) :: modname="tabfi"
[3455]90
[3184]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
[3455]100        time=0
[3184]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.3144621/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, &
[3477]134               "Missing value for periastron_day in def files!",1)
[3184]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        lmixmin=30
145        emin_turb=1.e-6
146! optical properties of polar caps and ground emissivity
147        emisice(:)=0
148        emissiv=0
149        iceradius(:)=1.e-6 ! mean scat radius of N2 snow
150        dtemisice(:)=0 !time scale for snow metamorphism
151        volcapa=1000000 ! volumetric heat capacity of subsurface
152
153      ELSE
154!-----------------------------------------------------------------------
155!  Initialization of physical constants by reading array tab_cntrl(:)
156!               which contains these parameters (nid != 0 case)
157!-----------------------------------------------------------------------
158! Read 'controle' array
159!
160
161       call get_var("controle",tab_cntrl,found)
162       if (.not.found) then
163         call abort_physic(modname,"Failed reading <controle> array",1)
164       else
165         write(*,*)'tabfi: tab_cntrl',tab_cntrl
166       endif
167!
168!  Initialization of some physical constants
169! informations on physics grid
170!      if(ngrid.ne.tab_cntrl(tab0+1)) then
171!         print*,'tabfi: WARNING !!! tab_cntrl(tab0+1).ne.ngrid'
172!         print*,tab_cntrl(tab0+1),ngrid
173!      endif
174      lmax = nint(tab_cntrl(tab0+2))
175      day_ini = tab_cntrl(tab0+3)
176      time = tab_cntrl(tab0+4)
177      write (*,*) 'IN tabfi day_ini=',day_ini
178! Informations about planet for dynamics and physics
179      rad = tab_cntrl(tab0+5)
180      omeg = tab_cntrl(tab0+6)
181      g = tab_cntrl(tab0+7)
182      mugaz = tab_cntrl(tab0+8)
183      rcp = tab_cntrl(tab0+9)
184      cpp=(8.314511/(mugaz/1000.0))/rcp
185      daysec = tab_cntrl(tab0+10)
186      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)
[3453]195! for paleoclimate
196      tpal = tab_cntrl(tab0+20)
197      adjust = tab_cntrl(tab0+21) ! for Triton albedo adjustment
198      ! lmixmin = tab_cntrl(tab0+20)
199      ! emin_turb = tab_cntrl(tab0+21)
200
201!!! AF24: parameters below are not used?
[3184]202! optical properties of polar caps and ground emissivity
203      emisice(1) = tab_cntrl(tab0+24)
204      emisice(2) = tab_cntrl(tab0+25)
205      emissiv    = tab_cntrl(tab0+26)
206      iceradius(1)= tab_cntrl(tab0+31) ! mean scat radius of N2 snow (north)
207      iceradius(2)= tab_cntrl(tab0+32) ! mean scat radius of N2 snow (south)
208      dtemisice(1)= tab_cntrl(tab0+33) !time scale for snow metamorphism (north)
209      dtemisice(2)= tab_cntrl(tab0+34) !time scale for snow metamorphism (south)
210! soil properties
211      volcapa = tab_cntrl(tab0+35) ! volumetric heat capacity
212!-----------------------------------------------------------------------
213!       Save some constants for later use (as routine arguments)
214!-----------------------------------------------------------------------
215      p_omeg = omeg
216      p_g = g
217      p_cpp = cpp
218      p_mugaz = mugaz
219      p_daysec = daysec
220      p_rad=rad
221
[3455]222      ENDIF    ! end of (nid = 0)
[3184]223
224!-----------------------------------------------------------------------
225!       Write physical constants to output before modifying them
226!-----------------------------------------------------------------------
[3455]227
[3184]228   6  FORMAT(a20,e15.6,e15.6)
229   5  FORMAT(a20,f12.2,f12.2)
[3455]230
[3184]231      write(*,*) '*****************************************************'
232      write(*,*) 'Reading tab_cntrl when calling tabfi before changes'
233      write(*,*) '*****************************************************'
234      write(*,5) '(1)        = ngrid?',tab_cntrl(tab0+1),float(ngrid)
235      write(*,5) '(2)            lmax',tab_cntrl(tab0+2),float(lmax)
236      write(*,5) '(3)         day_ini',tab_cntrl(tab0+3),float(day_ini)
237      write(*,5) '(5)             rad',tab_cntrl(tab0+5),rad
238      write(*,5) '(10)         daysec',tab_cntrl(tab0+10),daysec
239      write(*,6) '(6)            omeg',tab_cntrl(tab0+6),omeg
240      write(*,5) '(7)               g',tab_cntrl(tab0+7),g
241      write(*,5) '(8)           mugaz',tab_cntrl(tab0+8),mugaz
242      write(*,5) '(9)             rcp',tab_cntrl(tab0+9),rcp
243      write(*,6) '(11)        dtphys?',tab_cntrl(tab0+11),dtphys
244
245      write(*,5) '(14)       year_day',tab_cntrl(tab0+14),year_day
246      write(*,5) '(15)       periastr',tab_cntrl(tab0+15),periastr
247      write(*,5) '(16)        apoastr',tab_cntrl(tab0+16),apoastr
248      write(*,5) '(17)       peri_day',tab_cntrl(tab0+17),peri_day
249      write(*,5) '(18)       obliquit',tab_cntrl(tab0+18),obliquit
250
251      write(*,6) '(19)             z0',tab_cntrl(tab0+19),z0
252      write(*,6) '(21)      emin_turb',tab_cntrl(tab0+21),emin_turb
253      write(*,5) '(20)        lmixmin',tab_cntrl(tab0+20),lmixmin
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(*,*) '(21)       emin_turb :  minimal energy (PBL)'
288      write(*,*) '(20)         lmixmin : mixing length (PBL)'
289      write(*,*) '(26)         emissiv : ground emissivity'
290      write(*,*) '(24 et 25)   emisice : N2 ice max emissivity '
291      write(*,*) '(31 et 32) iceradius : mean scat radius of N2 snow'
292      write(*,*) '(33 et 34) dtemisice : time scale for snow metamorphism'
293      write(*,*) '(35)      volcapa : soil volumetric heat capacity'
294      write(*,*) '(18)     obliquit : planet obliquity (deg)'
295      write(*,*) '(17)     peri_day : periastron date (sols since Ls=0)'
296      write(*,*) '(15)     periastr : min. star-planet dist (UA)'
297      write(*,*) '(16)     apoastr  : max. star-planet (UA)'
298      write(*,*) '(14)     year_day : length of year (in sols)'
299      write(*,*) '(5)      rad      : radius of the planet (m)'
300      write(*,*) '(6)      omeg     : planet rotation rate (rad/s)'
301      write(*,*) '(7)      g        : gravity (m/s2)'
302      write(*,*) '(8)      mugaz    : molecular mass '
303      write(*,*) '                       of the atmosphere (g/mol)'
304      write(*,*) '(9)      rcp      : r/Cp'
305      write(*,*) '(8)+(9)  calc_cpp_mugaz : r/Cp and mugaz '
306      write(*,*) '                 computed from gases.def'
307      write(*,*) '(10)     daysec   : length of a sol (s)'
308      write(*,*)
[3455]309
310
[3184]311      do while(modif(1:1).ne.'hello')
312        write(*,*)
313        write(*,*)
314        write(*,*) 'Changes to perform ?'
315        write(*,*) '   (enter keyword or return )'
316        write(*,*)
317        read(*,fmt='(a20)') modif
318        if (modif(1:1) .eq. ' ') goto 999
[3455]319
[3184]320        write(*,*)
321        write(*,*) modif(1:len_trim(modif)) , ' : '
322
323        if (modif(1:len_trim(modif)) .eq. 'day_ini') then
324          write(*,*) 'current value:',day_ini
325          write(*,*) 'enter new value:'
326 101      read(*,*,iostat=ierr) day_ini
327          if(ierr.ne.0) goto 101
328          write(*,*) ' '
329          write(*,*) 'day_ini (new value):',day_ini
330
331        else if (modif(1:len_trim(modif)) .eq. 'z0') then
332          write(*,*) 'current value:',z0
333          write(*,*) 'enter new value:'
334 102      read(*,*,iostat=ierr) z0
335          if(ierr.ne.0) goto 102
336          write(*,*) ' '
337          write(*,*) ' z0 (new value):',z0
338
339        else if (modif(1:len_trim(modif)) .eq. 'emin_turb') then
340          write(*,*) 'current value:',emin_turb
341          write(*,*) 'enter new value:'
342 103      read(*,*,iostat=ierr) emin_turb
343          if(ierr.ne.0) goto 103
344          write(*,*) ' '
345          write(*,*) ' emin_turb (new value):',emin_turb
346
347        else if (modif(1:len_trim(modif)) .eq. 'lmixmin') then
348          write(*,*) 'current value:',lmixmin
349          write(*,*) 'enter new value:'
350 104      read(*,*,iostat=ierr) lmixmin
351          if(ierr.ne.0) goto 104
352          write(*,*) ' '
353          write(*,*) ' lmixmin (new value):',lmixmin
354
355        else if (modif(1:len_trim(modif)) .eq. 'emissiv') then
356          write(*,*) 'current value:',emissiv
357          write(*,*) 'enter new value:'
358 105      read(*,*,iostat=ierr) emissiv
359          if(ierr.ne.0) goto 105
360          write(*,*) ' '
361          write(*,*) ' emissiv (new value):',emissiv
362
363        else if (modif(1:len_trim(modif)) .eq. 'emisice') then
364          write(*,*) 'current value emisice(1) North:',emisice(1)
365          write(*,*) 'enter new value:'
366 106      read(*,*,iostat=ierr) emisice(1)
367          if(ierr.ne.0) goto 106
[3455]368          write(*,*)
[3184]369          write(*,*) ' emisice(1) (new value):',emisice(1)
370          write(*,*)
371
372          write(*,*) 'current value emisice(2) South:',emisice(2)
373          write(*,*) 'enter new value:'
374 107      read(*,*,iostat=ierr) emisice(2)
375          if(ierr.ne.0) goto 107
[3455]376          write(*,*)
[3184]377          write(*,*) ' emisice(2) (new value):',emisice(2)
378
379        else if (modif(1:len_trim(modif)) .eq. 'iceradius') then
380          write(*,*) 'current value iceradius(1) North:',iceradius(1)
381          write(*,*) 'enter new value:'
382 110      read(*,*,iostat=ierr) iceradius(1)
383          if(ierr.ne.0) goto 110
[3455]384          write(*,*)
[3184]385          write(*,*) ' iceradius(1) (new value):',iceradius(1)
386          write(*,*)
387
388          write(*,*) 'current value iceradius(2) South:',iceradius(2)
389          write(*,*) 'enter new value:'
390 111      read(*,*,iostat=ierr) iceradius(2)
391          if(ierr.ne.0) goto 111
[3455]392          write(*,*)
[3184]393          write(*,*) ' iceradius(2) (new value):',iceradius(2)
394
395        else if (modif(1:len_trim(modif)) .eq. 'dtemisice') then
396          write(*,*) 'current value dtemisice(1) North:',dtemisice(1)
397          write(*,*) 'enter new value:'
398 112      read(*,*,iostat=ierr) dtemisice(1)
399          if(ierr.ne.0) goto 112
[3455]400          write(*,*)
[3184]401          write(*,*) ' dtemisice(1) (new value):',dtemisice(1)
402          write(*,*)
403
404          write(*,*) 'current value dtemisice(2) South:',dtemisice(2)
405          write(*,*) 'enter new value:'
406 113      read(*,*,iostat=ierr) dtemisice(2)
407          if(ierr.ne.0) goto 113
[3455]408          write(*,*)
[3184]409          write(*,*) ' dtemisice(2) (new value):',dtemisice(2)
410
411        else if (modif(1:len_trim(modif)) .eq. 'obliquit') then
412          write(*,*) 'current value:',obliquit
413          write(*,*) 'obliquit should be 25.19 on current Mars'
414          write(*,*) 'enter new value:'
415 115      read(*,*,iostat=ierr) obliquit
416          if(ierr.ne.0) goto 115
[3455]417          write(*,*)
[3184]418          write(*,*) ' obliquit (new value):',obliquit
419
420        else if (modif(1:len_trim(modif)) .eq. 'peri_day') then
421          write(*,*) 'current value:',peri_day
422          write(*,*) 'peri_day should be 485 on current Mars'
423          write(*,*) 'enter new value:'
424 116      read(*,*,iostat=ierr) peri_day
425          if(ierr.ne.0) goto 116
[3455]426          write(*,*)
[3184]427          write(*,*) ' peri_day (new value):',peri_day
428
429        else if (modif(1:len_trim(modif)) .eq. 'periastr') then
430          write(*,*) 'current value:',periastr
431          write(*,*) 'periastr should be 1.3814 AU on present-day Mars'
432          write(*,*) 'enter new value:'
433 117      read(*,*,iostat=ierr) periastr
434          if(ierr.ne.0) goto 117
[3455]435          write(*,*)
[3184]436          write(*,*) ' periastr (new value):',periastr
[3455]437
[3184]438        else if (modif(1:len_trim(modif)) .eq. 'apoastr') then
439          write(*,*) 'current value:',apoastr
440          write(*,*) 'apoastr should be 1.666 AU on present-day Mars'
441          write(*,*) 'enter new value:'
442 118      read(*,*,iostat=ierr) apoastr
443          if(ierr.ne.0) goto 118
[3455]444          write(*,*)
[3184]445          write(*,*) ' apoastr (new value):',apoastr
[3455]446
[3184]447        else if (modif(1:len_trim(modif)) .eq. 'volcapa') then
448          write(*,*) 'current value:',volcapa
449          write(*,*) 'enter new value:'
450 119      read(*,*,iostat=ierr) volcapa
451          if(ierr.ne.0) goto 119
[3455]452          write(*,*)
[3184]453          write(*,*) ' volcapa (new value):',volcapa
[3455]454
[3184]455        else if (modif(1:len_trim(modif)).eq.'rad') then
456          write(*,*) 'current value:',rad
457          write(*,*) 'enter new value:'
458 120      read(*,*,iostat=ierr) rad
459          if(ierr.ne.0) goto 120
[3455]460          write(*,*)
[3184]461          write(*,*) ' rad (new value):',rad
462
463        else if (modif(1:len_trim(modif)).eq.'omeg') then
464          write(*,*) 'current value:',omeg
465          write(*,*) 'enter new value:'
466 121      read(*,*,iostat=ierr) omeg
467          if(ierr.ne.0) goto 121
[3455]468          write(*,*)
[3184]469          write(*,*) ' omeg (new value):',omeg
[3455]470
[3184]471        else if (modif(1:len_trim(modif)).eq.'g') then
472          write(*,*) 'current value:',g
473          write(*,*) 'enter new value:'
474 122      read(*,*,iostat=ierr) g
475          if(ierr.ne.0) goto 122
[3455]476          write(*,*)
[3184]477          write(*,*) ' g (new value):',g
478
479        else if (modif(1:len_trim(modif)).eq.'mugaz') then
480          write(*,*) 'current value:',mugaz
481          write(*,*) 'enter new value:'
482 123      read(*,*,iostat=ierr) mugaz
483          if(ierr.ne.0) goto 123
[3455]484          write(*,*)
[3184]485          write(*,*) ' mugaz (new value):',mugaz
486          r=8.314511/(mugaz/1000.0)
487          write(*,*) ' R (new value):',r
488
489        else if (modif(1:len_trim(modif)).eq.'rcp') then
490          write(*,*) 'current value:',rcp
491          write(*,*) 'enter new value:'
492 124      read(*,*,iostat=ierr) rcp
493          if(ierr.ne.0) goto 124
[3455]494          write(*,*)
[3184]495          write(*,*) ' rcp (new value):',rcp
496          r=8.314511/(mugaz/1000.0)
497          cpp=r/rcp
498          write(*,*) ' cpp (new value):',cpp
499
500        else if (modif(1:len_trim(modif)).eq.'calc_cpp_mugaz') then
501          write(*,*) 'current value rcp, mugaz:',rcp,mugaz
502          cpp_mugaz_mode = 2
503          call su_gases
504          call calc_cpp_mugaz
[3455]505          write(*,*)
[3184]506          write(*,*) ' cpp (new value):',cpp
507          write(*,*) ' mugaz (new value):',mugaz
508          r=8.314511/(mugaz/1000.0)
509          rcp=r/cpp
510          write(*,*) ' rcp (new value):',rcp
[3455]511
[3184]512        else if (modif(1:len_trim(modif)).eq.'daysec') then
513          write(*,*) 'current value:',daysec
514          write(*,*) 'enter new value:'
515 125      read(*,*,iostat=ierr) daysec
516          if(ierr.ne.0) goto 125
[3455]517          write(*,*)
[3184]518          write(*,*) ' daysec (new value):',daysec
519
520!         added by RW!
521        else if (modif(1:len_trim(modif)).eq.'year_day') then
522          write(*,*) 'current value:',year_day
[3455]523          write(*,*) 'enter new value:'
[3184]524 126      read(*,*,iostat=ierr) year_day
525          if(ierr.ne.0) goto 126
526          write(*,*)
527          write(*,*) ' year_day (new value):',year_day
528
529        endif
530      enddo ! of do while(modif(1:1).ne.'hello')
531
532 999  continue
533
534!----------------------------------------------------------------------
535!       Write values of physical constants after modifications
536!-----------------------------------------------------------------------
[3455]537
[3184]538      write(*,*) '*****************************************************'
539      write(*,*) 'Reading tab_cntrl when calling tabfi AFTER changes'
540      write(*,*) '*****************************************************'
541      write(*,5) '(1)        = ngrid?',tab_cntrl(tab0+1),float(ngrid)
542      write(*,5) '(2)            lmax',tab_cntrl(tab0+2),float(lmax)
543      write(*,5) '(3)         day_ini',tab_cntrl(tab0+3),float(day_ini)
544      write(*,5) '(5)             rad',tab_cntrl(tab0+5),rad
545      write(*,5) '(10)         daysec',tab_cntrl(tab0+10),daysec
546      write(*,6) '(6)            omeg',tab_cntrl(tab0+6),omeg
547      write(*,5) '(7)               g',tab_cntrl(tab0+7),g
548      write(*,5) '(8)           mugaz',tab_cntrl(tab0+8),mugaz
549      write(*,5) '(9)             rcp',tab_cntrl(tab0+9),rcp
550      write(*,6) '(11)        dtphys?',tab_cntrl(tab0+11),dtphys
[3455]551
[3184]552      write(*,5) '(14)       year_day',tab_cntrl(tab0+14),year_day
553      write(*,5) '(15)       periastr',tab_cntrl(tab0+15),periastr
554      write(*,5) '(16)        apoastr',tab_cntrl(tab0+16),apoastr
555      write(*,5) '(17)       peri_day',tab_cntrl(tab0+17),peri_day
556      write(*,5) '(18)       obliquit',tab_cntrl(tab0+18),obliquit
[3455]557
[3184]558      write(*,6) '(19)             z0',tab_cntrl(tab0+19),z0
559      write(*,6) '(21)      emin_turb',tab_cntrl(tab0+21),emin_turb
560      write(*,5) '(20)        lmixmin',tab_cntrl(tab0+20),lmixmin
[3455]561
[3184]562      write(*,5) '(26)        emissiv',tab_cntrl(tab0+26),emissiv
563      write(*,5) '(24)     emisice(1)',tab_cntrl(tab0+24),emisice(1)
564      write(*,5) '(25)     emisice(2)',tab_cntrl(tab0+25),emisice(2)
565      write(*,6) '(31)   iceradius(1)',tab_cntrl(tab0+31),iceradius(1)
566      write(*,6) '(32)   iceradius(2)',tab_cntrl(tab0+32),iceradius(2)
567      write(*,5) '(33)   dtemisice(1)',tab_cntrl(tab0+33),dtemisice(1)
568      write(*,5) '(34)   dtemisice(2)',tab_cntrl(tab0+34),dtemisice(2)
[3455]569
[3184]570      write(*,5) '(35)        volcapa',tab_cntrl(tab0+35),volcapa
571
[3455]572      write(*,*)
573      write(*,*)
[3184]574
575      ENDIF ! of if (Lmodif == 1)
576
577!-----------------------------------------------------------------------
578!       Save some constants for later use (as routine arguments)
579!-----------------------------------------------------------------------
580      p_omeg = omeg
581      p_g = g
582      p_cpp = cpp
583      p_mugaz = mugaz
584      p_daysec = daysec
585      p_rad=rad
586
587
588      END SUBROUTINE tabfi
589
590end module tabfi_mod
Note: See TracBrowser for help on using the repository browser.