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

Last change on this file since 310 was 310, checked in by millour, 10 years ago

Add the possibility to have a sponge acting on the temperature
in the physics, using flag "callradsponge=.true.". The sponge
vertical extention is the same as the sponge layer in the dynamics.
EM

File size: 28.3 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         if (is_master) write(*,*) 
207     &      "Add sponge after radiative transfer?"
208         callradsponge=.false.
209         call getin_p("callradsponge", callradsponge)
210         if (is_master) write(*,*) "callradsponge = ", callradsponge
211
212! Test of incompatibility:
213! if tlocked, then diurnal should be false
214         if (tlocked.and.diurnal) then
215           if (is_master) print*,'If diurnal=true, we should turn off',
216     &          ' tlocked.'
217           stop
218         endif
219
220         if (is_master) write(*,*) "Tidal resonance ratio ?"
221         nres=0          ! default value
222         call getin_p("nres",nres)
223         if (is_master) write(*,*) "nres = ",nres
224
225         if (is_master) write(*,*) "Write some extra output to the",
226     &          " screen ?"
227         lwrite=.false. ! default value
228         call getin_p("lwrite",lwrite)
229         if (is_master) write(*,*) " lwrite = ",lwrite
230
231         if (is_master) write(*,*) "Save statistics in file stats.nc ?"
232         callstats=.true. ! default value
233         call getin_p("callstats",callstats)
234         if (is_master) write(*,*) " callstats = ",callstats
235
236         if (is_master) write(*,*) "Test energy conservation of model",
237     &          " physics ?"
238         enertest=.false. ! default value
239         call getin_p("enertest",enertest)
240         if (is_master) write(*,*) " enertest = ",enertest
241
242         if (is_master) write(*,*) "Check to see if cpp values used",
243     &          " match gases.def ?"
244         check_cpp_match=.true. ! default value
245         call getin_p("check_cpp_match",check_cpp_match)
246         if (is_master) write(*,*) " check_cpp_match = ",check_cpp_match
247
248         if (is_master) write(*,*) "call radiative transfer ?"
249         callrad=.true. ! default value
250         call getin_p("callrad",callrad)
251         if (is_master) write(*,*) " callrad = ",callrad
252
253         if (is_master) write(*,*) "call correlated-k radiative",
254     &          " transfer ?"
255         corrk=.true. ! default value
256         call getin_p("corrk",corrk)
257         if (is_master) write(*,*) " corrk = ",corrk
258
259         if (is_master) write(*,*) "prohibit calculations outside",
260     &          " corrk T grid?"
261         strictboundcorrk=.true. ! default value
262         call getin_p("strictboundcorrk",strictboundcorrk)
263         if (is_master) write(*,*) "strictboundcorrk = ",
264     &          strictboundcorrk
265
266         if (is_master) write(*,*) "call gaseous absorption in the",
267     &          " visible bands? (matters only if callrad=T)"
268         callgasvis=.false. ! default value
269         call getin_p("callgasvis",callgasvis)
270         if (is_master) write(*,*) " callgasvis = ",callgasvis
271       
272         if (is_master) write(*,*) "call continuum opacities in",
273     &          " radiative transfer ? (matters only if callrad=T)"
274         continuum=.true. ! default value
275         call getin_p("continuum",continuum)
276         if (is_master) write(*,*) " continuum = ",continuum
277
278         if (is_master) write(*,*) "use analytic function for H2O",
279     &          " continuum ?"
280         H2Ocont_simple=.false. ! default value
281         call getin_p("H2Ocont_simple",H2Ocont_simple)
282         if (is_master) write(*,*) " H2Ocont_simple = ",H2Ocont_simple
283 
284         if (is_master) write(*,*) "call turbulent vertical",
285     &          " diffusion ?"
286         calldifv=.true. ! default value
287         call getin_p("calldifv",calldifv)
288         if (is_master) write(*,*) " calldifv = ",calldifv
289
290         if (is_master) write(*,*) "use turbdiff instead of vdifc ?"
291         UseTurbDiff=.true. ! default value
292         call getin_p("UseTurbDiff",UseTurbDiff)
293         if (is_master) write(*,*) " UseTurbDiff = ",UseTurbDiff
294
295         if (is_master) write(*,*) "call convective adjustment ?"
296         calladj=.true. ! default value
297         call getin_p("calladj",calladj)
298         if (is_master) write(*,*) " calladj = ",calladj
299
300         if (is_master) write(*,*) "call CO2 condensation ?"
301         co2cond=.false. ! default value
302         call getin_p("co2cond",co2cond)
303         if (is_master) write(*,*) " co2cond = ",co2cond
304! Test of incompatibility
305         if (co2cond.and.(.not.tracer)) then
306            if (is_master) print*,'We need a CO2 ice tracer to',
307     &          ' condense CO2'
308            call abort
309         endif
310 
311         if (is_master) write(*,*) "CO2 supersaturation level ?"
312         co2supsat=1.0 ! default value
313         call getin_p("co2supsat",co2supsat)
314         if (is_master) write(*,*) " co2supsat = ",co2supsat
315
316         if (is_master) write(*,*) "Radiative timescale for Newtonian",
317     &          " cooling ?"
318         tau_relax=30. ! default value
319         call getin_p("tau_relax",tau_relax)
320         if (is_master) write(*,*) " tau_relax = ",tau_relax
321         tau_relax=tau_relax*24*3600 ! convert Earth days --> seconds
322
323         if (is_master) write(*,*)"call thermal conduction in",
324     &          " the soil ?"
325         callsoil=.true. ! default value
326         call getin_p("callsoil",callsoil)
327         if (is_master) write(*,*) " callsoil = ",callsoil
328         
329         if (is_master) write(*,*)"Rad transfer is computed",
330     &          " every iradia physical timestep"
331         iradia=1 ! default value
332         call getin_p("iradia",iradia)
333         if (is_master) write(*,*)" iradia = ",iradia
334       
335         if (is_master) write(*,*)"Rayleigh scattering ?"
336         rayleigh=.false.
337         call getin_p("rayleigh",rayleigh)
338         if (is_master) write(*,*)" rayleigh = ",rayleigh
339
340         if (is_master) write(*,*) "Use blackbody for stellar",
341     &          " spectrum ?"
342         stelbbody=.false. ! default value
343         call getin_p("stelbbody",stelbbody)
344         if (is_master) write(*,*) " stelbbody = ",stelbbody
345
346         if (is_master) write(*,*) "Stellar blackbody temperature ?"
347         stelTbb=5800.0 ! default value
348         call getin_p("stelTbb",stelTbb)
349         if (is_master) write(*,*) " stelTbb = ",stelTbb
350
351         if (is_master) write(*,*)"Output mean OLR in 1D?"
352         meanOLR=.false.
353         call getin_p("meanOLR",meanOLR)
354         if (is_master) write(*,*)" meanOLR = ",meanOLR
355
356         if (is_master) write(*,*)"Output spectral OLR in 3D?"
357         specOLR=.false.
358         call getin_p("specOLR",specOLR)
359         if (is_master) write(*,*)" specOLR = ",specOLR
360
361         if (is_master) write(*,*)"Operate in kastprof mode?"
362         kastprof=.false.
363         call getin_p("kastprof",kastprof)
364         if (is_master) write(*,*)" kastprof = ",kastprof
365
366         if (is_master) write(*,*)"Uniform absorption in",
367     &          " radiative transfer?"
368         graybody=.false.
369         call getin_p("graybody",graybody)
370         if (is_master) write(*,*)" graybody = ",graybody
371
372! Slab Ocean
373         if (is_master) write(*,*) "Use slab-ocean ?"
374         ok_slab_ocean=.false.         ! default value
375         call getin_p("ok_slab_ocean",ok_slab_ocean)
376         if (is_master) write(*,*) "ok_slab_ocean = ",ok_slab_ocean
377
378         if (is_master) write(*,*) "Use slab-sea-ice ?"
379         ok_slab_sic=.true.         ! default value
380         call getin_p("ok_slab_sic",ok_slab_sic)
381         if (is_master) write(*,*) "ok_slab_sic = ",ok_slab_sic
382
383         if (is_master) write(*,*) "Use heat transport for the ocean ?"
384         ok_slab_heat_transp=.true.   ! default value
385         call getin_p("ok_slab_heat_transp",ok_slab_heat_transp)
386         if (is_master) write(*,*) "ok_slab_heat_transp = ",
387     &          ok_slab_heat_transp
388
389
390
391! Test of incompatibility:
392! if kastprof used, we must be in 1D
393         if (kastprof.and.(ngrid.gt.1)) then
394           if (is_master) print*,'kastprof can only be used in 1D!'
395           call abort
396         endif
397
398         if (is_master) write(*,*)"Stratospheric temperature",
399     &          " for kastprof mode?"
400         Tstrat=167.0
401         call getin_p("Tstrat",Tstrat)
402         if (is_master) write(*,*)" Tstrat = ",Tstrat
403
404         if (is_master) write(*,*)"Remove lower boundary?"
405         nosurf=.false.
406         call getin_p("nosurf",nosurf)
407         if (is_master) write(*,*)" nosurf = ",nosurf
408
409! Tests of incompatibility:
410         if (nosurf.and.callsoil) then
411           if (is_master) print*,'nosurf not compatible with soil',
412     &          ' scheme! ... got to make a choice!'
413           call abort
414         endif
415
416         if (is_master) write(*,*)"Add an internal heat flux?",
417     .             "... matters only if callsoil=F"
418         intheat=0.
419         call getin_p("intheat",intheat)
420         if (is_master) write(*,*)" intheat = ",intheat
421
422         if (is_master) write(*,*)"Use Newtonian cooling for",
423     &          " radiative transfer?"
424         newtonian=.false.
425         call getin_p("newtonian",newtonian)
426         if (is_master) write(*,*)" newtonian = ",newtonian
427
428! Tests of incompatibility:
429         if (newtonian.and.corrk) then
430           if (is_master) print*,'newtonian not compatible with ',
431     &          'correlated-k!'
432           call abort
433         endif
434         if (newtonian.and.calladj) then
435           if (is_master) print*,'newtonian not compatible with ',
436     &          'adjustment!'
437           call abort
438         endif
439         if (newtonian.and.calldifv) then
440           if (is_master) print*,'newtonian not compatible with a ',
441     &          'boundary layer!'
442           call abort
443         endif
444
445         if (is_master) write(*,*)"Test physics timescale in 1D?"
446         testradtimes=.false.
447         call getin_p("testradtimes",testradtimes)
448         if (is_master) write(*,*)" testradtimes = ",testradtimes
449
450! Test of incompatibility:
451! if testradtimes used, we must be in 1D
452         if (testradtimes.and.(ngrid.gt.1)) then
453           if (is_master) print*,'testradtimes can only be used in 1D!'
454           call abort
455         endif
456
457         if (is_master) write(*,*)"Default planetary temperature?"
458         tplanet=215.0
459         call getin_p("tplanet",tplanet)
460         if (is_master) write(*,*)" tplanet = ",tplanet
461
462         if (is_master) write(*,*)"Which star?"
463         startype=1 ! default value = Sol
464         call getin_p("startype",startype)
465         if (is_master) write(*,*)" startype = ",startype
466
467         if (is_master) write(*,*)"Value of stellar flux at 1 AU?"
468         Fat1AU=1356.0 ! default value = Sol today
469         call getin_p("Fat1AU",Fat1AU)
470         if (is_master) write(*,*)" Fat1AU = ",Fat1AU
471
472
473! TRACERS:
474
475         if (is_master) write(*,*)"Varying H2O cloud fraction?"
476         CLFvarying=.false.     ! default value
477         call getin_p("CLFvarying",CLFvarying)
478         if (is_master) write(*,*)" CLFvarying = ",CLFvarying
479
480         if (is_master) write(*,*)"Value of fixed H2O cloud fraction?"
481         CLFfixval=1.0                ! default value
482         call getin_p("CLFfixval",CLFfixval)
483         if (is_master) write(*,*)" CLFfixval = ",CLFfixval
484
485         if (is_master) write(*,*)"fixed radii for Cloud particles?"
486         radfixed=.false. ! default value
487         call getin_p("radfixed",radfixed)
488         if (is_master) write(*,*)" radfixed = ",radfixed
489
490         if(kastprof)then
491            radfixed=.true.
492         endif 
493
494         if (is_master) write(*,*)"Number mixing ratio of CO2 ice particles:"
495         Nmix_co2=1.e6 ! default value
496         call getin_p("Nmix_co2",Nmix_co2)
497         if (is_master) write(*,*)" Nmix_co2 = ",Nmix_co2
498
499!         write(*,*)"Number of radiatively active aerosols:"
500!         naerkind=0. ! default value
501!         call getin_p("naerkind",naerkind)
502!         write(*,*)" naerkind = ",naerkind
503
504         if (is_master) write(*,*)"Opacity of dust (if used):"
505         dusttau=0. ! default value
506         call getin_p("dusttau",dusttau)
507         if (is_master) write(*,*)" dusttau = ",dusttau
508
509         if (is_master) write(*,*)"Radiatively active CO2 aerosols?"
510         aeroco2=.false.     ! default value
511         call getin_p("aeroco2",aeroco2)
512         if (is_master) write(*,*)" aeroco2 = ",aeroco2
513
514         if (is_master) write(*,*)"Fixed CO2 aerosol distribution?"
515         aerofixco2=.false.     ! default value
516         call getin_p("aerofixco2",aerofixco2)
517         if (is_master) write(*,*)" aerofixco2 = ",aerofixco2
518
519         if (is_master) write(*,*)"Radiatively active water ice?"
520         aeroh2o=.false.     ! default value
521         call getin_p("aeroh2o",aeroh2o)
522         if (is_master) write(*,*)" aeroh2o = ",aeroh2o
523
524         if (is_master) write(*,*)"Fixed H2O aerosol distribution?"
525         aerofixh2o=.false.     ! default value
526         call getin_p("aerofixh2o",aerofixh2o)
527         if (is_master) write(*,*)" aerofixh2o = ",aerofixh2o
528
529         if (is_master) write(*,*)"Radiatively active sulfuric",
530     &          " acid aersols?"
531         aeroh2so4=.false.     ! default value
532         call getin_p("aeroh2so4",aeroh2so4)
533         if (is_master) write(*,*)" aeroh2so4 = ",aeroh2so4
534         
535!=================================
536
537         if (is_master) write(*,*)"Radiatively active two-layer aersols?"
538         aeroback2lay=.false.     ! default value
539         call getin_p("aeroback2lay",aeroback2lay)
540         if (is_master) write(*,*)" aeroback2lay = ",aeroback2lay
541
542         if (is_master) write(*,*)"TWOLAY AEROSOL: total optical depth ",
543     &              "in the tropospheric layer (visible)"
544         obs_tau_col_tropo=8.D0
545         call getin_p("obs_tau_col_tropo",obs_tau_col_tropo)
546         if (is_master) write(*,*)" obs_tau_col_tropo = ",
547     &          obs_tau_col_tropo
548
549         if (is_master) write(*,*)"TWOLAY AEROSOL: total optical depth ",
550     &              "in the stratospheric layer (visible)"
551         obs_tau_col_strato=0.08D0
552         call getin_p("obs_tau_col_strato",obs_tau_col_strato)
553         if (is_master) write(*,*)" obs_tau_col_strato = ",
554     &          obs_tau_col_strato
555
556         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_bottom_tropo?",
557     &          " in pa"
558         pres_bottom_tropo=66000.0
559         call getin_p("pres_bottom_tropo",pres_bottom_tropo)
560         if (is_master) write(*,*)" pres_bottom_tropo = ",
561     &          pres_bottom_tropo
562
563         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_top_tropo?",
564     &          " in pa"
565         pres_top_tropo=18000.0
566         call getin_p("pres_top_tropo",pres_top_tropo)
567         if (is_master) write(*,*)" pres_top_tropo = ",pres_top_tropo
568
569         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_bottom_strato?",
570     &          " in pa"
571         pres_bottom_strato=2000.0
572         call getin_p("pres_bottom_strato",pres_bottom_strato)
573         if (is_master) write(*,*)" pres_bottom_strato = ",
574     &          pres_bottom_strato
575
576         if (is_master) write(*,*)"TWOLAY AEROSOL: pres_top_strato?",
577     &          " in pa"
578         pres_top_strato=100.0
579         call getin_p("pres_top_strato",pres_top_strato)
580         if (is_master) write(*,*)" pres_top_strato = ",pres_top_strato
581
582         if (is_master) write(*,*)"TWOLAY AEROSOL: particle size",
583     &          " in the tropospheric layer, in meters"
584         size_tropo=2.e-6
585         call getin_p("size_tropo",size_tropo)
586         if (is_master) write(*,*)" size_tropo = ",size_tropo
587
588         if (is_master) write(*,*)"TWOLAY AEROSOL: particle size",
589     &           "in the stratospheric layer, in meters"
590         size_strato=1.e-7
591         call getin_p("size_strato",size_strato)
592         if (is_master) write(*,*)" size_strato = ",size_strato
593
594!=================================
595
596         if (is_master) write(*,*)"Cloud pressure level (with ",
597     &          "kastprof only):"
598         cloudlvl=0. ! default value
599         call getin_p("cloudlvl",cloudlvl)
600         if (is_master) write(*,*)" cloudlvl = ",cloudlvl
601
602         if (is_master) write(*,*)"Is the variable gas species",
603     &          " radiatively active?"
604         Tstrat=167.0
605         varactive=.false.
606         call getin_p("varactive",varactive)
607         if (is_master) write(*,*)" varactive = ",varactive
608
609         if (is_master) write(*,*)"Is the variable gas species",
610     &          " distribution set?"
611         varfixed=.false.
612         call getin_p("varfixed",varfixed)
613         if (is_master) write(*,*)" varfixed = ",varfixed
614
615         if (is_master) write(*,*)"What is the saturation % of",
616     &          " the variable species?"
617         satval=0.8
618         call getin_p("satval",satval)
619         if (is_master) write(*,*)" satval = ",satval
620
621
622! Test of incompatibility:
623! if varactive, then varfixed should be false
624         if (varactive.and.varfixed) then
625           if (is_master) print*,'if varactive, varfixed must be OFF!'
626           stop
627         endif
628
629         if (is_master) write(*,*) "Gravitationnal sedimentation ?"
630         sedimentation=.false. ! default value
631         call getin_p("sedimentation",sedimentation)
632         if (is_master) write(*,*) " sedimentation = ",sedimentation
633
634         if (is_master) write(*,*) "Compute water cycle ?"
635         water=.false. ! default value
636         call getin_p("water",water)
637         if (is_master) write(*,*) " water = ",water
638         
639! Test of incompatibility:
640! if water is true, there should be at least a tracer
641         if (water.and.(.not.tracer)) then
642           if (is_master) print*,'if water is ON, tracer must be ',
643     &          'ON too!'
644           stop
645         endif
646
647         if (is_master) write(*,*) "Include water condensation ?"
648         watercond=.false. ! default value
649         call getin_p("watercond",watercond)
650         if (is_master) write(*,*) " watercond = ",watercond
651
652! Test of incompatibility:
653! if watercond is used, then water should be used too
654         if (watercond.and.(.not.water)) then
655           if (is_master) print*,'if watercond is used, water should ',
656     &          'be used too'
657           stop
658         endif
659
660         if (is_master) write(*,*) "Include water precipitation ?"
661         waterrain=.false. ! default value
662         call getin_p("waterrain",waterrain)
663         if (is_master) write(*,*) " waterrain = ",waterrain
664
665         if (is_master) write(*,*) "Include surface hydrology ?"
666         hydrology=.false. ! default value
667         call getin_p("hydrology",hydrology)
668         if (is_master) write(*,*) " hydrology = ",hydrology
669
670         if (is_master) write(*,*) "Evolve surface water sources ?"
671         sourceevol=.false. ! default value
672         call getin_p("sourceevol",sourceevol)
673         if (is_master) write(*,*) " sourceevol = ",sourceevol
674
675         if (is_master) write(*,*) "Ice evolution timestep ?"
676         icetstep=100.0 ! default value
677         call getin_p("icetstep",icetstep)
678         if (is_master) write(*,*) " icetstep = ",icetstep
679
680         if (is_master) write(*,*) "Snow albedo ?"
681         albedosnow=0.5         ! default value
682         call getin_p("albedosnow",albedosnow)
683         if (is_master) write(*,*) " albedosnow = ",albedosnow
684
685         if (is_master) write(*,*) "Maximum ice thickness ?"
686         maxicethick=2.0         ! default value
687         call getin_p("maxicethick",maxicethick)
688         if (is_master) write(*,*) " maxicethick = ",maxicethick
689
690         if (is_master) write(*,*) "Freezing point of seawater ?"
691         Tsaldiff=-1.8          ! default value
692         call getin_p("Tsaldiff",Tsaldiff)
693         if (is_master) write(*,*) " Tsaldiff = ",Tsaldiff
694
695         if (is_master) write(*,*) "Does user want to force",
696     &          " cpp and mugaz?"
697         force_cpp=.false. ! default value
698         call getin_p("force_cpp",force_cpp)
699         if (is_master) write(*,*) " force_cpp = ",force_cpp
700
701         if (force_cpp) then
702           mugaz = -99999.
703           if (is_master) PRINT *,'MEAN MOLECULAR MASS in g mol-1 ?'
704           call getin_p("mugaz",mugaz)
705           IF (mugaz.eq.-99999.) THEN
706               if (is_master) PRINT *, "mugaz must be set if ',
707     &          'force_cpp = T"
708               STOP
709           ELSE
710               if (is_master) write(*,*) "inifis: mugaz=",mugaz
711           ENDIF
712           !Ehouarn: once mugaz has been set, r, the specific
713           ! gas constant must be computed
714           r=8.3144622/(mugaz*1.e-3)
715           if (is_master) write(*,*) "inifis: r=",r
716           
717           cpp = -99999.
718           if (is_master) PRINT *,'SPECIFIC HEAT CAPACITY in ',
719     &          'J K-1 kg-1 ?'
720           call getin_p("cpp",cpp)
721           IF (cpp.eq.-99999.) THEN
722               if (is_master) PRINT *, "cpp must be set if ',
723     &          'force_cpp = T"
724               STOP
725           ELSE
726               if (is_master) write(*,*) "inifis: cpp=",cpp
727           ENDIF
728!         else
729!           mugaz=8.314*1000./pr
730         endif
731         call su_gases
732         call calc_cpp_mugaz
733
734         if (is_master) PRINT*,'------------------------------------'
735         if (is_master) PRINT*
736         if (is_master) PRINT*
737      ELSE
738         if (is_master) write(*,*)
739         if (is_master) write(*,*) 'Cannot read file callphys.def.",
740     &          " Is it here ?'
741         stop
742      ENDIF
743
7448000  FORMAT(t5,a12,l8)
7458001  FORMAT(t5,a12,i8)
746
747      if (is_master) then
748      PRINT*
749      PRINT*,'inifis: daysec',daysec
750      PRINT*
751      PRINT*,'inifis: The radiative transfer is computed:'
752      PRINT*,'           each ',iradia,' physical time-step'
753      PRINT*,'        or each ',iradia*dtphys,' seconds'
754      PRINT*
755      end if
756
757
758!-----------------------------------------------------------------------
759!     Some more initialization:
760!     ------------------------
761
762      ! ALLOCATE ARRAYS IN comgeomfi_h
763      IF (.not. ALLOCATED(lati)) ALLOCATE(lati(ngrid))
764      IF (.not. ALLOCATED(long)) ALLOCATE(long(ngrid))
765      IF (.not. ALLOCATED(area)) ALLOCATE(area(ngrid))
766
767      CALL SCOPY(ngrid,plon,1,long,1)
768      CALL SCOPY(ngrid,plat,1,lati,1)
769      CALL SCOPY(ngrid,parea,1,area,1)
770      totarea=SSUM(ngrid,area,1)
771      call planetwide_sumval(area,totarea_planet)
772
773      !! those are defined in comdiurn_h.F90
774      IF (.not.ALLOCATED(sinlat)) ALLOCATE(sinlat(ngrid))
775      IF (.not.ALLOCATED(coslat)) ALLOCATE(coslat(ngrid))
776      IF (.not.ALLOCATED(sinlon)) ALLOCATE(sinlon(ngrid))
777      IF (.not.ALLOCATED(coslon)) ALLOCATE(coslon(ngrid))
778
779      DO ig=1,ngrid
780         sinlat(ig)=sin(plat(ig))
781         coslat(ig)=cos(plat(ig))
782         sinlon(ig)=sin(plon(ig))
783         coslon(ig)=cos(plon(ig))
784      ENDDO
785
786!$OMP MASTER
787      pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h
788!$OMP END MASTER
789!$OMP BARRIER
790
791      ! allocate "comsoil_h" arrays
792      call ini_comsoil_h(ngrid)
793     
794      END
Note: See TracBrowser for help on using the repository browser.