source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/inifis.F @ 298

Last change on this file since 298 was 298, checked in by milmd, 10 years ago

Less output messages are written. On 20000 cores it is better. In LMDZ, only master of MPI and OpenMP can write.

File size: 28.0 KB
Line 
1      SUBROUTINE inifis(ngrid,nlayer,nq,
2     $           day_ini,pdaysec,ptimestep,
3     $           plat,plon,parea,
4     $           prad,pg,pr,pcpp)
5
6      use radinc_h, only : naerkind
7      use datafile_mod, only: datadir
8      use comdiurn_h, only: sinlat, coslat, sinlon, coslon
9      use comgeomfi_h, only: long, lati, area, totarea, totarea_planet
10      use comsoil_h, only: ini_comsoil_h
11      use control_mod, only: ecritphy
12      use planete_mod, only: nres
13      use planetwide_mod, only: planetwide_sumval
14
15!=======================================================================
16!
17!   purpose:
18!   -------
19!
20!   Initialisation for the physical parametrisations of the LMD
21!   Generic Model.
22!
23!   author: Frederic Hourdin 15 / 10 /93
24!   -------
25!   modified: Sebastien Lebonnois 11/06/2003 (new callphys.def)
26!             Ehouarn Millour (oct. 2008) tracers are now identified
27!              by their names and may not be contiguously
28!              stored in the q(:,:,:,:) array
29!             E.M. (june 2009) use getin routine to load parameters
30!
31!
32!   arguments:
33!   ----------
34!
35!   input:
36!   ------
37!
38!    ngrid                 Size of the horizontal grid.
39!                          All internal loops are performed on that grid.
40!    nlayer                Number of vertical layers.
41!    pdayref               Day of reference for the simulation
42!    pday                  Number of days counted from the North. Spring
43!                          equinoxe.
44!
45!=======================================================================
46!
47!-----------------------------------------------------------------------
48!   declarations:
49!   -------------
50      use datafile_mod, only: datadir
51! to use  'getin'
52!      USE ioipsl_getincom
53      USE ioipsl_getincom_p
54      use mod_phys_lmdz_para, only : is_master
55      IMPLICIT NONE
56!#include "dimensions.h"
57!#include "dimphys.h"
58!#include "planete.h"
59#include "comcstfi.h"
60#include "callkeys.h"
61
62
63
64      REAL,INTENT(IN) :: prad,pg,pr,pcpp,pdaysec,ptimestep
65 
66      INTEGER,INTENT(IN) :: ngrid,nlayer,nq
67      REAL,INTENT(IN) :: plat(ngrid),plon(ngrid),parea(ngrid)
68      integer day_ini
69      INTEGER ig,ierr
70 
71      EXTERNAL iniorbit,orbite
72      EXTERNAL SSUM
73      REAL SSUM
74 
75      CHARACTER ch1*12
76      CHARACTER ch80*80
77
78      logical chem, h2o
79      logical :: parameter, doubleq=.false.
80
81      real psurf,pN2 ! added by RW for Gliese 581d N2+CO2
82
83!$OMP MASTER
84      rad=prad
85      daysec=pdaysec
86      dtphys=ptimestep
87      cpp=pcpp
88      g=pg
89      r=pr
90      rcp=r/cpp
91      ! Ehouarn debug:
92      if (is_master) then
93      write(*,*) "inifis: rad=",rad
94      write(*,*) "        daysec=",daysec
95      write(*,*) "        dtphys=",dtphys
96      write(*,*) "        cpp=",cpp
97      write(*,*) "        g=",g
98      write(*,*) "        r=",r
99      write(*,*) "        rcp=",rcp
100      end if
101
102      avocado = 6.02214179e23   ! added by RW
103
104!$OMP END MASTER
105!$OMP BARRIER
106
107      ! read in 'ecritphy' (frequency of calls to physics, in dynamical steps)
108      ! (also done in dyn3d/defrun_new but not in LMDZ.COMMON)
109      call getin_p("ecritphy",ecritphy)
110
111! --------------------------------------------------------------
112!  Reading the "callphys.def" file controlling some key options
113! --------------------------------------------------------------
114     
115!$OMP MASTER     
116      ! check that 'callphys.def' file is around
117      OPEN(99,file='callphys.def',status='old',form='formatted'
118     &     ,iostat=ierr)
119      CLOSE(99)
120      IF(ierr.EQ.0) iscallphys=.true. !iscallphys initialised as false in callkeys.h
121!$OMP END MASTER
122!$OMP BARRIER
123     
124!!!      IF(ierr.EQ.0) THEN
125      IF(iscallphys) THEN
126         if (is_master) then
127         PRINT*
128         PRINT*
129         PRINT*,'--------------------------------------------'
130         PRINT*,' inifis: Parametres pour la physique (callphys.def)'
131         PRINT*,'--------------------------------------------'
132         end if
133
134         if (is_master) write(*,*) "Directory where external input",
135     &          " files are:"
136
137         ! default 'datadir' is set in "datadir_mod"
138         call getin_p("datadir",datadir) ! default path
139         if (is_master) write(*,*) " datadir = ",trim(datadir)
140
141         if (is_master) write(*,*) "Run with or without tracer",
142     &          " transport ?"
143         tracer=.false. ! default value
144         call getin_p("tracer",tracer)
145         if (is_master) write(*,*) " tracer = ",tracer
146
147         if (is_master) write(*,*) "Run with or without atm mass",
148     &          " update due to tracer evaporation/condensation?"
149         mass_redistrib=.false. ! default value
150         call getin_p("mass_redistrib",mass_redistrib)
151         if (is_master) write(*,*) " mass_redistrib = ",mass_redistrib
152
153         if (is_master) write(*,*) "Diurnal cycle ?"
154         if (is_master) write(*,*) "(if diurnal=false, diurnal",
155     &          " averaged solar heating)"
156         diurnal=.true. ! default value
157         call getin_p("diurnal",diurnal)
158         if (is_master) write(*,*) " diurnal = ",diurnal
159
160         if (is_master) write(*,*) "Seasonal cycle ?"
161         if (is_master) write(*,*) "(if season=false, Ls stays",
162     &          " constant, to value set in 'start'"
163         season=.true. ! default value
164         call getin_p("season",season)
165         if (is_master) write(*,*) " season = ",season
166
167         if (is_master) write(*,*) "Tidally resonant rotation ?"
168         tlocked=.false. ! default value
169         call getin_p("tlocked",tlocked)
170         if (is_master) write(*,*) "tlocked = ",tlocked
171
172         if (is_master) write(*,*) "Saturn ring shadowing ?"
173         rings_shadow = .false.
174         call getin_p("rings_shadow", rings_shadow)
175         if (is_master) write(*,*) "rings_shadow = ", rings_shadow
176         
177         if (is_master) write(*,*) "Compute latitude-dependent",
178     &          " gravity field?"
179         oblate = .false.
180         call getin_p("oblate", oblate)
181         if (is_master) write(*,*) "oblate = ", oblate
182
183         if (is_master) write(*,*) "Flattening of the planet (a-b)/a "
184         flatten = 0.0
185         call getin_p("flatten", flatten)
186         if (is_master) write(*,*) "flatten = ", flatten
187         
188
189         if (is_master) write(*,*) "Needed if oblate=.true.: J2"
190         J2 = 0.0
191         call getin_p("J2", J2)
192         if (is_master) write(*,*) "J2 = ", J2
193         
194         if (is_master) write(*,*) "Needed if oblate=.true.: Planet",
195     &          " mass (*1e24 kg)"
196         MassPlanet = 0.0
197         call getin_p("MassPlanet", MassPlanet)
198         if (is_master) write(*,*) "MassPlanet = ", MassPlanet         
199
200         if (is_master) write(*,*) "Needed if oblate=.true.: Planet",
201     &          " mean radius (m)"
202         Rmean = 0.0
203         call getin_p("Rmean", Rmean)
204         if (is_master) write(*,*) "Rmean = ", Rmean
205         
206! Test of incompatibility:
207! if tlocked, then diurnal should be false
208         if (tlocked.and.diurnal) then
209           if (is_master) print*,'If diurnal=true, we should turn off',
210     &          ' tlocked.'
211           stop
212         endif
213
214         if (is_master) write(*,*) "Tidal resonance ratio ?"
215         nres=0          ! default value
216         call getin_p("nres",nres)
217         if (is_master) write(*,*) "nres = ",nres
218
219         if (is_master) write(*,*) "Write some extra output to the",
220     &          " screen ?"
221         lwrite=.false. ! default value
222         call getin_p("lwrite",lwrite)
223         if (is_master) write(*,*) " lwrite = ",lwrite
224
225         if (is_master) write(*,*) "Save statistics in file stats.nc ?"
226         callstats=.true. ! default value
227         call getin_p("callstats",callstats)
228         if (is_master) write(*,*) " callstats = ",callstats
229
230         if (is_master) write(*,*) "Test energy conservation of model",
231     &          " physics ?"
232         enertest=.false. ! default value
233         call getin_p("enertest",enertest)
234         if (is_master) write(*,*) " enertest = ",enertest
235
236         if (is_master) write(*,*) "Check to see if cpp values used",
237     &          " match gases.def ?"
238         check_cpp_match=.true. ! default value
239         call getin_p("check_cpp_match",check_cpp_match)
240         if (is_master) write(*,*) " check_cpp_match = ",check_cpp_match
241
242         if (is_master) write(*,*) "call radiative transfer ?"
243         callrad=.true. ! default value
244         call getin_p("callrad",callrad)
245         if (is_master) write(*,*) " callrad = ",callrad
246
247         if (is_master) write(*,*) "call correlated-k radiative",
248     &          " transfer ?"
249         corrk=.true. ! default value
250         call getin_p("corrk",corrk)
251         if (is_master) write(*,*) " corrk = ",corrk
252
253         if (is_master) write(*,*) "prohibit calculations outside",
254     &          " corrk T grid?"
255         strictboundcorrk=.true. ! default value
256         call getin_p("strictboundcorrk",strictboundcorrk)
257         if (is_master) write(*,*) "strictboundcorrk = ",
258     &          strictboundcorrk
259
260         if (is_master) write(*,*) "call gaseous absorption in the",
261     &          " visible bands? (matters only if callrad=T)"
262         callgasvis=.false. ! default value
263         call getin_p("callgasvis",callgasvis)
264         if (is_master) write(*,*) " callgasvis = ",callgasvis
265       
266         if (is_master) write(*,*) "call continuum opacities in",
267     &          " radiative transfer ? (matters only if callrad=T)"
268         continuum=.true. ! default value
269         call getin_p("continuum",continuum)
270         if (is_master) write(*,*) " continuum = ",continuum
271
272         if (is_master) write(*,*) "use analytic function for H2O",
273     &          " continuum ?"
274         H2Ocont_simple=.false. ! default value
275         call getin_p("H2Ocont_simple",H2Ocont_simple)
276         if (is_master) write(*,*) " H2Ocont_simple = ",H2Ocont_simple
277 
278         if (is_master) write(*,*) "call turbulent vertical",
279     &          " diffusion ?"
280         calldifv=.true. ! default value
281         call getin_p("calldifv",calldifv)
282         if (is_master) write(*,*) " calldifv = ",calldifv
283
284         if (is_master) write(*,*) "use turbdiff instead of vdifc ?"
285         UseTurbDiff=.true. ! default value
286         call getin_p("UseTurbDiff",UseTurbDiff)
287         if (is_master) write(*,*) " UseTurbDiff = ",UseTurbDiff
288
289         if (is_master) write(*,*) "call convective adjustment ?"
290         calladj=.true. ! default value
291         call getin_p("calladj",calladj)
292         if (is_master) write(*,*) " calladj = ",calladj
293
294         if (is_master) write(*,*) "call CO2 condensation ?"
295         co2cond=.false. ! default value
296         call getin_p("co2cond",co2cond)
297         if (is_master) write(*,*) " co2cond = ",co2cond
298! Test of incompatibility
299         if (co2cond.and.(.not.tracer)) then
300            if (is_master) print*,'We need a CO2 ice tracer to',
301     &          ' condense CO2'
302            call abort
303         endif
304 
305         if (is_master) write(*,*) "CO2 supersaturation level ?"
306         co2supsat=1.0 ! default value
307         call getin_p("co2supsat",co2supsat)
308         if (is_master) write(*,*) " co2supsat = ",co2supsat
309
310         if (is_master) write(*,*) "Radiative timescale for Newtonian",
311     &          " cooling ?"
312         tau_relax=30. ! default value
313         call getin_p("tau_relax",tau_relax)
314         if (is_master) write(*,*) " tau_relax = ",tau_relax
315         tau_relax=tau_relax*24*3600 ! convert Earth days --> seconds
316
317         if (is_master) write(*,*)"call thermal conduction in",
318     &          " the soil ?"
319         callsoil=.true. ! default value
320         call getin_p("callsoil",callsoil)
321         if (is_master) write(*,*) " callsoil = ",callsoil
322         
323         if (is_master) write(*,*)"Rad transfer is computed",
324     &          " every iradia physical timestep"
325         iradia=1 ! default value
326         call getin_p("iradia",iradia)
327         if (is_master) write(*,*)" iradia = ",iradia
328       
329         if (is_master) write(*,*)"Rayleigh scattering ?"
330         rayleigh=.false.
331         call getin_p("rayleigh",rayleigh)
332         if (is_master) write(*,*)" rayleigh = ",rayleigh
333
334         if (is_master) write(*,*) "Use blackbody for stellar",
335     &          " spectrum ?"
336         stelbbody=.false. ! default value
337         call getin_p("stelbbody",stelbbody)
338         if (is_master) write(*,*) " stelbbody = ",stelbbody
339
340         if (is_master) write(*,*) "Stellar blackbody temperature ?"
341         stelTbb=5800.0 ! default value
342         call getin_p("stelTbb",stelTbb)
343         if (is_master) write(*,*) " stelTbb = ",stelTbb
344
345         if (is_master) write(*,*)"Output mean OLR in 1D?"
346         meanOLR=.false.
347         call getin_p("meanOLR",meanOLR)
348         if (is_master) write(*,*)" meanOLR = ",meanOLR
349
350         if (is_master) write(*,*)"Output spectral OLR in 3D?"
351         specOLR=.false.
352         call getin_p("specOLR",specOLR)
353         if (is_master) write(*,*)" specOLR = ",specOLR
354
355         if (is_master) write(*,*)"Operate in kastprof mode?"
356         kastprof=.false.
357         call getin_p("kastprof",kastprof)
358         if (is_master) write(*,*)" kastprof = ",kastprof
359
360         if (is_master) write(*,*)"Uniform absorption in",
361     &          " radiative transfer?"
362         graybody=.false.
363         call getin_p("graybody",graybody)
364         if (is_master) write(*,*)" graybody = ",graybody
365
366! Slab Ocean
367         if (is_master) write(*,*) "Use slab-ocean ?"
368         ok_slab_ocean=.false.         ! default value
369         call getin_p("ok_slab_ocean",ok_slab_ocean)
370         if (is_master) write(*,*) "ok_slab_ocean = ",ok_slab_ocean
371
372         if (is_master) write(*,*) "Use slab-sea-ice ?"
373         ok_slab_sic=.true.         ! default value
374         call getin_p("ok_slab_sic",ok_slab_sic)
375         if (is_master) write(*,*) "ok_slab_sic = ",ok_slab_sic
376
377         if (is_master) write(*,*) "Use heat transport for the ocean ?"
378         ok_slab_heat_transp=.true.   ! default value
379         call getin_p("ok_slab_heat_transp",ok_slab_heat_transp)
380         if (is_master) write(*,*) "ok_slab_heat_transp = ",
381     &          ok_slab_heat_transp
382
383
384
385! Test of incompatibility:
386! if kastprof used, we must be in 1D
387         if (kastprof.and.(ngrid.gt.1)) then
388           if (is_master) print*,'kastprof can only be used in 1D!'
389           call abort
390         endif
391
392         if (is_master) write(*,*)"Stratospheric temperature",
393     &          " for kastprof mode?"
394         Tstrat=167.0
395         call getin_p("Tstrat",Tstrat)
396         if (is_master) write(*,*)" Tstrat = ",Tstrat
397
398         if (is_master) write(*,*)"Remove lower boundary?"
399         nosurf=.false.
400         call getin_p("nosurf",nosurf)
401         if (is_master) write(*,*)" nosurf = ",nosurf
402
403! Tests of incompatibility:
404         if (nosurf.and.callsoil) then
405           if (is_master) print*,'nosurf not compatible with soil',
406     &          ' scheme! ... got to make a choice!'
407           call abort
408         endif
409
410         if (is_master) write(*,*)"Add an internal heat flux?",
411     .             "... matters only if callsoil=F"
412         intheat=0.
413         call getin_p("intheat",intheat)
414         if (is_master) write(*,*)" intheat = ",intheat
415
416         if (is_master) write(*,*)"Use Newtonian cooling for",
417     &          " radiative transfer?"
418         newtonian=.false.
419         call getin_p("newtonian",newtonian)
420         if (is_master) write(*,*)" newtonian = ",newtonian
421
422! Tests of incompatibility:
423         if (newtonian.and.corrk) then
424           if (is_master) print*,'newtonian not compatible with ',
425     &          'correlated-k!'
426           call abort
427         endif
428         if (newtonian.and.calladj) then
429           if (is_master) print*,'newtonian not compatible with ',
430     &          'adjustment!'
431           call abort
432         endif
433         if (newtonian.and.calldifv) then
434           if (is_master) print*,'newtonian not compatible with a ',
435     &          'boundary layer!'
436           call abort
437         endif
438
439         if (is_master) write(*,*)"Test physics timescale in 1D?"
440         testradtimes=.false.
441         call getin_p("testradtimes",testradtimes)
442         if (is_master) write(*,*)" testradtimes = ",testradtimes
443
444! Test of incompatibility:
445! if testradtimes used, we must be in 1D
446         if (testradtimes.and.(ngrid.gt.1)) then
447           if (is_master) print*,'testradtimes can only be used in 1D!'
448           call abort
449         endif
450
451         if (is_master) write(*,*)"Default planetary temperature?"
452         tplanet=215.0
453         call getin_p("tplanet",tplanet)
454         if (is_master) write(*,*)" tplanet = ",tplanet
455
456         if (is_master) write(*,*)"Which star?"
457         startype=1 ! default value = Sol
458         call getin_p("startype",startype)
459         if (is_master) write(*,*)" startype = ",startype
460
461         if (is_master) write(*,*)"Value of stellar flux at 1 AU?"
462         Fat1AU=1356.0 ! default value = Sol today
463         call getin_p("Fat1AU",Fat1AU)
464         if (is_master) write(*,*)" Fat1AU = ",Fat1AU
465
466
467! TRACERS:
468
469         if (is_master) write(*,*)"Varying H2O cloud fraction?"
470         CLFvarying=.false.     ! default value
471         call getin_p("CLFvarying",CLFvarying)
472         if (is_master) write(*,*)" CLFvarying = ",CLFvarying
473
474         if (is_master) write(*,*)"Value of fixed H2O cloud fraction?"
475         CLFfixval=1.0                ! default value
476         call getin_p("CLFfixval",CLFfixval)
477         if (is_master) write(*,*)" CLFfixval = ",CLFfixval
478
479         if (is_master) write(*,*)"fixed radii for Cloud particles?"
480         radfixed=.false. ! default value
481         call getin_p("radfixed",radfixed)
482         if (is_master) write(*,*)" radfixed = ",radfixed
483
484         if(kastprof)then
485            radfixed=.true.
486         endif 
487
488         if (is_master) write(*,*)"Number mixing ratio of CO2 ice particles:"
489         Nmix_co2=1.e6 ! default value
490         call getin_p("Nmix_co2",Nmix_co2)
491         if (is_master) write(*,*)" Nmix_co2 = ",Nmix_co2
492
493!         write(*,*)"Number of radiatively active aerosols:"
494!         naerkind=0. ! default value
495!         call getin_p("naerkind",naerkind)
496!         write(*,*)" naerkind = ",naerkind
497
498         if (is_master) write(*,*)"Opacity of dust (if used):"
499         dusttau=0. ! default value
500         call getin_p("dusttau",dusttau)
501         if (is_master) write(*,*)" dusttau = ",dusttau
502
503         if (is_master) write(*,*)"Radiatively active CO2 aerosols?"
504         aeroco2=.false.     ! default value
505         call getin_p("aeroco2",aeroco2)
506         if (is_master) write(*,*)" aeroco2 = ",aeroco2
507
508         if (is_master) write(*,*)"Fixed CO2 aerosol distribution?"
509         aerofixco2=.false.     ! default value
510         call getin_p("aerofixco2",aerofixco2)
511         if (is_master) write(*,*)" aerofixco2 = ",aerofixco2
512
513         if (is_master) write(*,*)"Radiatively active water ice?"
514         aeroh2o=.false.     ! default value
515         call getin_p("aeroh2o",aeroh2o)
516         if (is_master) write(*,*)" aeroh2o = ",aeroh2o
517
518         if (is_master) write(*,*)"Fixed H2O aerosol distribution?"
519         aerofixh2o=.false.     ! default value
520         call getin_p("aerofixh2o",aerofixh2o)
521         if (is_master) write(*,*)" aerofixh2o = ",aerofixh2o
522
523         if (is_master) write(*,*)"Radiatively active sulfuric",
524     &          " acid aersols?"
525         aeroh2so4=.false.     ! default value
526         call getin_p("aeroh2so4",aeroh2so4)
527         if (is_master) write(*,*)" aeroh2so4 = ",aeroh2so4
528         
529!=================================
530
531         if (is_master) write(*,*)"Radiatively active two-layer aersols?"
532         aeroback2lay=.false.     ! default value
533         call getin_p("aeroback2lay",aeroback2lay)
534         if (is_master) write(*,*)" aeroback2lay = ",aeroback2lay
535
536         if (is_master) write(*,*)"TWOLAY AEROSOL: total optical depth ",
537     &              "in the tropospheric layer (visible)"
538         obs_tau_col_tropo=8.D0
539         call getin_p("obs_tau_col_tropo",obs_tau_col_tropo)
540         if (is_master) write(*,*)" obs_tau_col_tropo = ",
541     &          obs_tau_col_tropo
542
543         if (is_master) write(*,*)"TWOLAY AEROSOL: total optical depth ",
544     &              "in the stratospheric layer (visible)"
545         obs_tau_col_strato=0.08D0
546         call getin_p("obs_tau_col_strato",obs_tau_col_strato)
547         if (is_master) write(*,*)" obs_tau_col_strato = ",
548     &          obs_tau_col_strato
549
550         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_bottom_tropo?",
551     &          " in pa"
552         pres_bottom_tropo=66000.0
553         call getin_p("pres_bottom_tropo",pres_bottom_tropo)
554         if (is_master) write(*,*)" pres_bottom_tropo = ",
555     &          pres_bottom_tropo
556
557         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_top_tropo?",
558     &          " in pa"
559         pres_top_tropo=18000.0
560         call getin_p("pres_top_tropo",pres_top_tropo)
561         if (is_master) write(*,*)" pres_top_tropo = ",pres_top_tropo
562
563         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_bottom_strato?",
564     &          " in pa"
565         pres_bottom_strato=2000.0
566         call getin_p("pres_bottom_strato",pres_bottom_strato)
567         if (is_master) write(*,*)" pres_bottom_strato = ",
568     &          pres_bottom_strato
569
570         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_top_strato?",
571     &          " in pa"
572         pres_top_strato=100.0
573         call getin_p("pres_top_strato",pres_top_strato)
574         if (is_master) write(*,*)" pres_top_strato = ",pres_top_strato
575
576         if (is_master) write(*,*)"TWOLAY AEROSOL: particle size",
577     &          " in the tropospheric layer, in meters"
578         size_tropo=2.e-6
579         call getin_p("size_tropo",size_tropo)
580         if (is_master) write(*,*)" size_tropo = ",size_tropo
581
582         if (is_master) write(*,*)"TWOLAY AEROSOL: particle size",
583     &           "in the stratospheric layer, in meters"
584         size_strato=1.e-7
585         call getin_p("size_strato",size_strato)
586         if (is_master) write(*,*)" size_strato = ",size_strato
587
588!=================================
589
590         if (is_master) write(*,*)"Cloud pressure level (with ",
591     &          "kastprof only):"
592         cloudlvl=0. ! default value
593         call getin_p("cloudlvl",cloudlvl)
594         if (is_master) write(*,*)" cloudlvl = ",cloudlvl
595
596         if (is_master) write(*,*)"Is the variable gas species",
597     &          " radiatively active?"
598         Tstrat=167.0
599         varactive=.false.
600         call getin_p("varactive",varactive)
601         if (is_master) write(*,*)" varactive = ",varactive
602
603         if (is_master) write(*,*)"Is the variable gas species",
604     &          " distribution set?"
605         varfixed=.false.
606         call getin_p("varfixed",varfixed)
607         if (is_master) write(*,*)" varfixed = ",varfixed
608
609         if (is_master) write(*,*)"What is the saturation % of",
610     &          " the variable species?"
611         satval=0.8
612         call getin_p("satval",satval)
613         if (is_master) write(*,*)" satval = ",satval
614
615
616! Test of incompatibility:
617! if varactive, then varfixed should be false
618         if (varactive.and.varfixed) then
619           if (is_master) print*,'if varactive, varfixed must be OFF!'
620           stop
621         endif
622
623         if (is_master) write(*,*) "Gravitationnal sedimentation ?"
624         sedimentation=.false. ! default value
625         call getin_p("sedimentation",sedimentation)
626         if (is_master) write(*,*) " sedimentation = ",sedimentation
627
628         if (is_master) write(*,*) "Compute water cycle ?"
629         water=.false. ! default value
630         call getin_p("water",water)
631         if (is_master) write(*,*) " water = ",water
632         
633! Test of incompatibility:
634! if water is true, there should be at least a tracer
635         if (water.and.(.not.tracer)) then
636           if (is_master) print*,'if water is ON, tracer must be ',
637     &          'ON too!'
638           stop
639         endif
640
641         if (is_master) write(*,*) "Include water condensation ?"
642         watercond=.false. ! default value
643         call getin_p("watercond",watercond)
644         if (is_master) write(*,*) " watercond = ",watercond
645
646! Test of incompatibility:
647! if watercond is used, then water should be used too
648         if (watercond.and.(.not.water)) then
649           if (is_master) print*,'if watercond is used, water should ',
650     &          'be used too'
651           stop
652         endif
653
654         if (is_master) write(*,*) "Include water precipitation ?"
655         waterrain=.false. ! default value
656         call getin_p("waterrain",waterrain)
657         if (is_master) write(*,*) " waterrain = ",waterrain
658
659         if (is_master) write(*,*) "Include surface hydrology ?"
660         hydrology=.false. ! default value
661         call getin_p("hydrology",hydrology)
662         if (is_master) write(*,*) " hydrology = ",hydrology
663
664         if (is_master) write(*,*) "Evolve surface water sources ?"
665         sourceevol=.false. ! default value
666         call getin_p("sourceevol",sourceevol)
667         if (is_master) write(*,*) " sourceevol = ",sourceevol
668
669         if (is_master) write(*,*) "Ice evolution timestep ?"
670         icetstep=100.0 ! default value
671         call getin_p("icetstep",icetstep)
672         if (is_master) write(*,*) " icetstep = ",icetstep
673
674         if (is_master) write(*,*) "Snow albedo ?"
675         albedosnow=0.5         ! default value
676         call getin_p("albedosnow",albedosnow)
677         if (is_master) write(*,*) " albedosnow = ",albedosnow
678
679         if (is_master) write(*,*) "Maximum ice thickness ?"
680         maxicethick=2.0         ! default value
681         call getin_p("maxicethick",maxicethick)
682         if (is_master) write(*,*) " maxicethick = ",maxicethick
683
684         if (is_master) write(*,*) "Freezing point of seawater ?"
685         Tsaldiff=-1.8          ! default value
686         call getin_p("Tsaldiff",Tsaldiff)
687         if (is_master) write(*,*) " Tsaldiff = ",Tsaldiff
688
689         if (is_master) write(*,*) "Does user want to force",
690     &          " cpp and mugaz?"
691         force_cpp=.false. ! default value
692         call getin_p("force_cpp",force_cpp)
693         if (is_master) write(*,*) " force_cpp = ",force_cpp
694
695         if (force_cpp) then
696           mugaz = -99999.
697           if (is_master) PRINT *,'MEAN MOLECULAR MASS in g mol-1 ?'
698           call getin_p("mugaz",mugaz)
699           IF (mugaz.eq.-99999.) THEN
700               if (is_master) PRINT *, "mugaz must be set if ',
701     &          'force_cpp = T"
702               STOP
703           ELSE
704               if (is_master) write(*,*) "inifis: mugaz=",mugaz
705           ENDIF
706           !Ehouarn: once mugaz has been set, r, the specific
707           ! gas constant must be computed
708           r=8.3144622/(mugaz*1.e-3)
709           if (is_master) write(*,*) "inifis: r=",r
710           
711           cpp = -99999.
712           if (is_master) PRINT *,'SPECIFIC HEAT CAPACITY in ',
713     &          'J K-1 kg-1 ?'
714           call getin_p("cpp",cpp)
715           IF (cpp.eq.-99999.) THEN
716               if (is_master) PRINT *, "cpp must be set if ',
717     &          'force_cpp = T"
718               STOP
719           ELSE
720               if (is_master) write(*,*) "inifis: cpp=",cpp
721           ENDIF
722!         else
723!           mugaz=8.314*1000./pr
724         endif
725         call su_gases
726         call calc_cpp_mugaz
727
728         if (is_master) PRINT*,'------------------------------------'
729         if (is_master) PRINT*
730         if (is_master) PRINT*
731      ELSE
732         if (is_master) write(*,*)
733         if (is_master) write(*,*) 'Cannot read file callphys.def.",
734     &          " Is it here ?'
735         stop
736      ENDIF
737
7388000  FORMAT(t5,a12,l8)
7398001  FORMAT(t5,a12,i8)
740
741      if (is_master) then
742      PRINT*
743      PRINT*,'inifis: daysec',daysec
744      PRINT*
745      PRINT*,'inifis: The radiative transfer is computed:'
746      PRINT*,'           each ',iradia,' physical time-step'
747      PRINT*,'        or each ',iradia*dtphys,' seconds'
748      PRINT*
749      end if
750
751
752!-----------------------------------------------------------------------
753!     Some more initialization:
754!     ------------------------
755
756      ! ALLOCATE ARRAYS IN comgeomfi_h
757      IF (.not. ALLOCATED(lati)) ALLOCATE(lati(ngrid))
758      IF (.not. ALLOCATED(long)) ALLOCATE(long(ngrid))
759      IF (.not. ALLOCATED(area)) ALLOCATE(area(ngrid))
760
761      CALL SCOPY(ngrid,plon,1,long,1)
762      CALL SCOPY(ngrid,plat,1,lati,1)
763      CALL SCOPY(ngrid,parea,1,area,1)
764      totarea=SSUM(ngrid,area,1)
765      call planetwide_sumval(area,totarea_planet)
766
767      !! those are defined in comdiurn_h.F90
768      IF (.not.ALLOCATED(sinlat)) ALLOCATE(sinlat(ngrid))
769      IF (.not.ALLOCATED(coslat)) ALLOCATE(coslat(ngrid))
770      IF (.not.ALLOCATED(sinlon)) ALLOCATE(sinlon(ngrid))
771      IF (.not.ALLOCATED(coslon)) ALLOCATE(coslon(ngrid))
772
773      DO ig=1,ngrid
774         sinlat(ig)=sin(plat(ig))
775         coslat(ig)=cos(plat(ig))
776         sinlon(ig)=sin(plon(ig))
777         coslon(ig)=cos(plon(ig))
778      ENDDO
779
780!$OMP MASTER
781      pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h
782!$OMP END MASTER
783!$OMP BARRIER
784
785      ! allocate "comsoil_h" arrays
786      call ini_comsoil_h(ngrid)
787     
788      END
Note: See TracBrowser for help on using the repository browser.