source: trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90 @ 3809

Last change on this file since 3809 was 3725, checked in by emillour, 3 months ago

Generic PCM:
Change the way the rate of outputs for diagfi.nc files is specified:
IMPORTANT: Specifying "ecritphy" no longer possible and will trigger an error.
Use "diagfi_output_rate" to specify output rate (in physics steps) instead.
This should makes things (hopefully) clearer for users and also better
enforces a cleaner and clearer separation between dynamics and physics.
EM

File size: 51.2 KB
Line 
1MODULE inifis_mod
2IMPLICIT NONE
3
4CONTAINS
5
6  SUBROUTINE inifis(ngrid,nlayer,nq, &
7             day_ini,pdaysec,nday,ptimestep, &
8             plat,plon,parea, &
9             prad,pg,pr,pcpp)
10
11  use init_print_control_mod, only: init_print_control
12  use radinc_h, only: ini_radinc_h, naerkind
13  use radcommon_h, only: ini_radcommon_h
14  use radii_mod, only: radfixed, Nmix_co2
15  use datafile_mod, only: datadir
16  use comdiurn_h, only: sinlat, coslat, sinlon, coslon
17  use comgeomfi_h, only: totarea, totarea_planet
18  use comsoil_h, only: ini_comsoil_h, nsoilmx, lay1_soil, alpha_soil
19  use time_phylmdz_mod, only: diagfi_output_rate, &
20                              init_time, daysec, dtphys
21  use comcstfi_mod, only: rad, cpp, g, r, rcp, &
22                          mugaz, pi, avocado
23  use planete_mod, only: nres
24  use planetwide_mod, only: planetwide_sumval
25  use callkeys_mod
26  use wstats_mod, only: callstats
27  use ioipsl_getin_p_mod, only : getin_p
28  use mod_phys_lmdz_para, only : is_parallel, is_master, bcast
29  use newton_cooling_hotJ, only: planetary_suffix
30
31!=======================================================================
32!
33!   purpose:
34!   -------
35!
36!   Initialisation for the physical parametrisations of the LMD
37!   Generic Model.
38!
39!   author: Frederic Hourdin 15 / 10 /93
40!   -------
41!   modified: Sebastien Lebonnois 11/06/2003 (new callphys.def)
42!             Ehouarn Millour (oct. 2008) tracers are now identified
43!              by their names and may not be contiguously
44!              stored in the q(:,:,:,:) array
45!             E.M. (june 2009) use getin routine to load parameters
46!
47!
48!   arguments:
49!   ----------
50!
51!   input:
52!   ------
53!
54!    ngrid                 Size of the horizontal grid.
55!                          All internal loops are performed on that grid.
56!    nlayer                Number of vertical layers.
57!    pdayref               Day of reference for the simulation
58!    pday                  Number of days counted from the North. Spring
59!                          equinoxe.
60!
61!=======================================================================
62!
63!-----------------------------------------------------------------------
64!   declarations:
65!   -------------
66  use datafile_mod, only: datadir
67  use ioipsl_getin_p_mod, only: getin_p
68  IMPLICIT NONE
69
70
71
72  REAL,INTENT(IN) :: prad,pg,pr,pcpp,pdaysec,ptimestep
73  INTEGER,INTENT(IN) :: nday
74  INTEGER,INTENT(IN) :: ngrid,nlayer,nq
75  REAL,INTENT(IN) :: plat(ngrid),plon(ngrid),parea(ngrid)
76  integer,intent(in) :: day_ini
77
78  INTEGER :: ig
79  REAL :: ecritphy ! to check that this obsolete flag is no longer used...
80  CHARACTER(len=20) :: rname="inifis" ! routine name, for messages
81 
82  EXTERNAL iniorbit,orbite
83  EXTERNAL SSUM
84  REAL SSUM
85 
86  ! Initialize flags lunout, prt_level, debug (in print_control_mod)
87  CALL init_print_control
88
89  ! initialize constants in comcstfi_mod
90  rad=prad
91  cpp=pcpp
92  g=pg
93  r=pr
94  rcp=r/cpp
95  mugaz=8.314463*1000./pr ! dummy init
96  pi=2.*asin(1.)
97  avocado = 6.022141e23   ! added by RW
98
99  ! Initialize some "temporal and calendar" related variables
100#ifndef MESOSCALE
101  CALL init_time(day_ini,pdaysec,nday,ptimestep)
102#endif
103
104  ! read in some parameters from "run.def" for physics,
105  ! or shared between dynamics and physics.
106  ecritphy=-666 ! dummy default value
107  call getin_p("ecritphy",ecritphy) ! frequency of outputs in physics,
108                                    ! in dynamical steps
109  if (ecritphy/=-666) then
110    call abort_physic(rname, &
111         "Error: parameter ecritphy is obsolete! Remove it! "//&
112         "And use diagfi_output_rate instead",1)
113  endif
114
115  ! do we read a startphy.nc file? (default: .true.)
116  call getin_p("startphy_file",startphy_file)
117 
118! --------------------------------------------------------------
119!  Reading the "callphys.def" file controlling some key options
120! --------------------------------------------------------------
121     
122  IF (is_master) THEN
123    ! check if 'callphys.def' file is around
124    inquire(file="callphys.def",exist=iscallphys)
125    write(*,*) trim(rname)//": iscallphys=",iscallphys
126  ENDIF
127  call bcast(iscallphys)
128
129  IF(iscallphys) THEN
130
131     if (is_master) then
132       write(*,*)
133       write(*,*)'--------------------------------------------'
134       write(*,*)trim(rname)//': Parameters for the physics (callphys.def)'
135       write(*,*)'--------------------------------------------'
136     endif
137
138     if (is_master) write(*,*) trim(rname)//&
139       ": Directory where external input files are:"
140     ! default 'datadir' is set in "datadir_mod"
141     call getin_p("datadir",datadir) ! default path
142     if (is_master) write(*,*) trim(rname)//": datadir = ",trim(datadir)
143
144     if (is_master) write(*,*) trim(rname)//&
145       ": Run with or without tracer transport ?"
146     tracer=.false. ! default value
147     call getin_p("tracer",tracer)
148     if (is_master) write(*,*) trim(rname)//": tracer = ",tracer
149
150     if (is_master) write(*,*) trim(rname)//&
151       ": Output rate for diagfi.nc file (in physics steps) ?"
152     diagfi_output_rate=24 !default value
153     call getin_p("diagfi_output_rate",diagfi_output_rate)
154     if (is_master) write(*,*) trim(rname)//": diagfi_output_rate = ",&
155                               diagfi_output_rate
156
157     if (is_master) write(*,*) trim(rname)//&
158       ": Run with or without atm mass update "//&
159       " due to tracer evaporation/condensation?"
160     mass_redistrib=.false. ! default value
161     call getin_p("mass_redistrib",mass_redistrib)
162     if (is_master) write(*,*) trim(rname)//&
163       ": mass_redistrib = ",mass_redistrib
164
165     if (is_master) then
166       write(*,*) trim(rname)//": Diurnal cycle ?"
167       write(*,*) trim(rname)//&
168       ": (if diurnal=false, diurnal averaged solar heating)"
169     endif
170     diurnal=.true. ! default value
171     call getin_p("diurnal",diurnal)
172     if (is_master) write(*,*) trim(rname)//": diurnal = ",diurnal
173
174     if (is_master) then
175       write(*,*) trim(rname)//": Seasonal cycle ?"
176       write(*,*) trim(rname)//&
177       ": (if season=false, Ls stays constant, to value ", &
178       "set in 'start'"
179     endif
180     season=.true. ! default value
181     call getin_p("season",season)
182     if (is_master) write(*,*) trim(rname)//": season = ",season
183     
184     if (is_master) write(*,*) trim(rname)//&
185       ": No seasonal cycle: initial day to lock the run during restart"
186     noseason_day=0.0 ! default value
187     call getin_p("noseason_day",noseason_day)
188     if (is_master) write(*,*) trim(rname)//": noseason_day=", noseason_day
189
190     if (is_master) write(*,*) trim(rname)//": Tidally resonant rotation ?"
191     tlocked=.false. ! default value
192     call getin_p("tlocked",tlocked)
193     if (is_master) write(*,*) trim(rname)//": tlocked = ",tlocked
194
195     if (is_master) write(*,*) trim(rname)//": Saturn ring shadowing ?"
196     rings_shadow = .false.
197     call getin_p("rings_shadow", rings_shadow)
198     if (is_master) write(*,*) trim(rname)//": rings_shadow = ", rings_shadow
199         
200     if (is_master) write(*,*) trim(rname)//&
201       ": Compute latitude-dependent gravity field?"
202     oblate = .false.
203     call getin_p("oblate", oblate)
204     if (is_master) write(*,*) trim(rname)//": oblate = ", oblate
205
206     if (is_master) write(*,*) trim(rname)//": Flattening of the planet (a-b)/a "
207     flatten = 0.0
208     call getin_p("flatten", flatten)
209     if (is_master) write(*,*) trim(rname)//": flatten = ", flatten
210         
211
212     if (is_master) write(*,*) trim(rname)//": Needed if oblate=.true.: J2"
213     J2 = 0.0
214     call getin_p("J2", J2)
215     if (is_master) write(*,*) trim(rname)//": J2 = ", J2
216         
217     if (is_master) write(*,*) trim(rname)//&
218       ": Needed if oblate=.true.: Planet mass (*1e24 kg)"
219     MassPlanet = 0.0
220     call getin_p("MassPlanet", MassPlanet)
221     if (is_master) write(*,*) trim(rname)//": MassPlanet = ", MassPlanet         
222
223     if (is_master) write(*,*) trim(rname)//&
224       ": Needed if oblate=.true.: Planet mean radius (m)"
225     Rmean = 0.0
226     call getin_p("Rmean", Rmean)
227     if (is_master) write(*,*) trim(rname)//": Rmean = ", Rmean
228         
229! Test of incompatibility:
230! if tlocked, then diurnal should be false
231     if (tlocked.and.diurnal) then
232       call abort_physic(rname,"If diurnal=true, we should turn off tlocked",1)
233     endif
234
235     if (is_master) write(*,*) trim(rname)//": Tidal resonance ratio ?"
236     nres=0          ! default value
237     call getin_p("nres",nres)
238     if (is_master) write(*,*) trim(rname)//": nres = ",nres
239
240     if (is_master) write(*,*) trim(rname)//&
241       ": Write some extra output to the screen ?"
242     lwrite=.false. ! default value
243     call getin_p("lwrite",lwrite)
244     if (is_master) write(*,*) trim(rname)//": lwrite = ",lwrite
245
246     if (is_master) write(*,*) trim(rname)//&
247       ": Save statistics in file stats.nc ?"
248     callstats=.true. ! default value
249     call getin_p("callstats",callstats)
250     if (is_master) write(*,*) trim(rname)//": callstats = ",callstats
251
252     if (is_master) write(*,*) trim(rname)//&
253       ": Test energy conservation of model physics ?"
254     enertest=.false. ! default value
255     call getin_p("enertest",enertest)
256     if (is_master) write(*,*) trim(rname)//": enertest = ",enertest
257
258     if (is_master) write(*,*) trim(rname)//": call radiative transfer ?"
259     callrad=.true. ! default value
260     call getin_p("callrad",callrad)
261     if (is_master) write(*,*) trim(rname)//": callrad = ",callrad
262
263     if (is_master) write(*,*) trim(rname)//&
264       ": call correlated-k radiative transfer ?"
265     corrk=.true. ! default value
266     call getin_p("corrk",corrk)
267     if (is_master) write(*,*) trim(rname)//": corrk = ",corrk
268
269     if (is_master) write(*,*) trim(rname)//&
270       ": prohibit calculations outside corrk T grid?"
271     strictboundcorrk=.true. ! default value
272     call getin_p("strictboundcorrk",strictboundcorrk)
273     if (is_master) write(*,*) trim(rname)//&
274       ": strictboundcorrk = ",strictboundcorrk
275       
276     if (is_master) write(*,*) trim(rname)//&
277       ": call thermosphere ?"
278     callthermos=.false. ! default value
279     call getin_p("callthermos",callthermos)
280     if (is_master) write(*,*) trim(rname)//&
281       ": callthermos = ",callthermos
282     if(callthermos) then
283       if (is_master) write(*,*) trim(rname)//&
284         ": flux from thermosphere ? W/m2"
285       phitop = 0.0 ! default value
286       call getin_p("phitop",phitop)
287       if (is_master) write(*,*) trim(rname)//&
288         ": phitop = ",phitop
289       if (is_master) write(*,*) trim(rname)//&
290         ": Thickness of conduction ? m"
291       zztop = 1000.
292       call getin_p("zztop",zztop)
293       if (is_master) write(*,*) trim(rname)//&
294         ": zztop = ",zztop
295       if (is_master) write(*,*) trim(rname)//&
296         ": Force conduction ? "
297       force_conduction=.false. ! default value
298       call getin_p("force_conduction",force_conduction)
299       if (is_master) write(*,*) trim(rname)//&
300         ": force_conduction = ",force_conduction
301       if(force_conduction) then
302         if (is_master) write(*,*) trim(rname)//&
303           ": a_coeff ?"
304         a_coeff = 0.0 ! default value
305         call getin_p("a_coeff",a_coeff)
306         if (is_master) write(*,*) trim(rname)//&
307           ": a_coeff = ",a_coeff
308         if (is_master) write(*,*) trim(rname)//&
309           ": s_coeff ?"
310         s_coeff = 0.0 ! default value
311         call getin_p("s_coeff",s_coeff)
312         if (is_master) write(*,*) trim(rname)//&
313           ": s_coeff = ",s_coeff
314       endif
315     endif
316
317     if (is_master) write(*,*) trim(rname)//&
318       ": prohibit calculations outside CIA T grid?"
319     strictboundcia=.true. ! default value
320     call getin_p("strictboundcia",strictboundcia)
321     if (is_master) write(*,*) trim(rname)//&
322       ": strictboundcia = ",strictboundcia
323
324     if (is_master) write(*,*) trim(rname)//&
325       ": Minimum atmospheric temperature for Planck function integration ?"
326     tplanckmin=30.0 ! default value
327     call getin_p("tplanckmin",tplanckmin)
328     if (is_master) write(*,*) trim(rname)//": tplanckmin = ",tplanckmin
329 
330     if (is_master) write(*,*) trim(rname)//&
331       ": Maximum atmospheric temperature for Planck function integration ?"
332     tplanckmax=1500.0 ! default value
333     call getin_p("tplanckmax",tplanckmax)
334     if (is_master) write(*,*) trim(rname)//": tplanckmax = ",tplanckmax
335 
336     if (is_master) write(*,*) trim(rname)//&
337       ": Temperature step for Planck function integration ?"
338     dtplanck=0.1 ! default value
339     call getin_p("dtplanck",dtplanck)
340     if (is_master) write(*,*) trim(rname)//": dtplanck = ",dtplanck
341 
342     if (is_master) write(*,*) trim(rname)//&
343       ": call gaseous absorption in the visible bands?"//&
344       " (matters only if callrad=T)"
345     callgasvis=.false. ! default value
346     call getin_p("callgasvis",callgasvis)
347     if (is_master) write(*,*) trim(rname)//": callgasvis = ",callgasvis
348       
349     if (is_master) write(*,*) trim(rname)//&
350       ": call continuum opacities in radiative transfer ?"//&
351       " (matters only if callrad=T)"
352     continuum=.true. ! default value
353     call getin_p("continuum",continuum)
354     if (is_master) write(*,*) trim(rname)//": continuum = ",continuum
355     
356     if (is_master) write(*,*) trim(rname)//&
357       ": use generic continuum opacity database ?"//&
358       " (matters only if callrad=T)"
359     generic_continuum_database=.true. ! default value
360     call getin_p("generic_continuum_database",generic_continuum_database)
361     if (is_master) write(*,*) trim(rname)//": generic_continuum_database = ", &
362     generic_continuum_database
363 
364     if (is_master) write(*,*) trim(rname)//": version for H2H2 CIA file ?"
365     versH2H2cia=2011 ! default value (should be 2018 but retrocompatibility first)
366     call getin_p("versH2H2cia",versH2H2cia)
367     if (is_master) write(*,*) trim(rname)//": versH2H2cia = ",versH2H2cia
368     ! Sanity check
369     if (versH2H2cia.ne.2011 .and. versH2H2cia.ne.2018) then
370        call abort_physic(rname,"Error: Choose a correct value (2011 or 2018) for versH2H2cia !",1)
371     endif
372     
373     if (is_master) write(*,*) trim(rname)//&
374       ": CIA - normal or equilibrium ortho-para mixture for H2?"
375     H2orthopara_mixture = 'normal'
376     call getin_p("H2orthopara_mixture",H2orthopara_mixture)
377     if (is_master) write(*,*)trim(rname)//&
378       ": H2orthopara_mixture = ",trim(H2orthopara_mixture)
379
380     if (is_master) write(*,*) trim(rname)//&
381       ": call turbulent vertical diffusion ?"
382     calldifv=.true. ! default value
383     call getin_p("calldifv",calldifv)
384     if (is_master) write(*,*) trim(rname)//": calldifv = ",calldifv
385
386     if (is_master) write(*,*) trim(rname)//": use turbdiff instead of vdifc ?"
387     UseTurbDiff=.true. ! default value
388     call getin_p("UseTurbDiff",UseTurbDiff)
389     if (is_master) write(*,*) trim(rname)//": UseTurbDiff = ",UseTurbDiff
390
391     if (is_master) write(*,*) trim(rname)//": call convective adjustment ?"
392     calladj=.true. ! default value
393     call getin_p("calladj",calladj)
394     if (is_master) write(*,*) trim(rname)//": calladj = ",calladj
395
396     if (is_master) write(*,*)trim(rname)//&
397       ": call Lott's non-oro GWs parameterisation scheme ?"
398     calllott_nonoro=.false. ! default value
399     call getin_p("calllott_nonoro",calllott_nonoro)
400     if (is_master) write(*,*)trim(rname)//&
401       ": calllott_nonoro = ",calllott_nonoro
402     
403     if (is_master) write(*,*) trim(rname)//": call thermal plume model ?"
404     calltherm=.false. ! default value
405     call getin_p("calltherm",calltherm)
406     if (is_master) write(*,*) trim(rname)//": calltherm = ",calltherm
407
408     if (is_master) write(*,*) trim(rname)//": call CO2 condensation ?"
409     co2cond=.false. ! default value
410     call getin_p("co2cond",co2cond)
411     if (is_master) write(*,*) trim(rname)//": co2cond = ",co2cond
412! Test of incompatibility
413     if (co2cond.and.(.not.tracer)) then
414        call abort_physic(rname,"Error: We need a CO2 ice tracer to condense CO2!",1)
415     endif
416 
417     if (is_master) write(*,*) trim(rname)//": CO2 supersaturation level ?"
418     co2supsat=1.0 ! default value
419     call getin_p("co2supsat",co2supsat)
420     if (is_master) write(*,*) trim(rname)//": co2supsat = ",co2supsat
421
422     if (is_master) write(*,*) trim(rname)//&
423     ": Radiative timescale for Newtonian cooling (in Earth days)?"
424     tau_relax=30. ! default value
425     call getin_p("tau_relax",tau_relax)
426     if (is_master) write(*,*) trim(rname)//": tau_relax = ",tau_relax
427     tau_relax=tau_relax*24*3600 ! convert Earth days --> seconds
428
429     if (is_master) write(*,*) trim(rname)//&
430     ": planetary suffix for loading file in Newtonian cooling HotJ ?"
431     planetary_suffix='default_value' ! default value
432     call getin_p("planetary_suffix",planetary_suffix)
433     if (is_master) write(*,*) trim(rname)//": planetary_suffix",planetary_suffix
434     
435     if (is_master) write(*,*)trim(rname)//&
436       ": call thermal conduction in the soil ?"
437     callsoil=.true. ! default value
438     call getin_p("callsoil",callsoil)
439     if (is_master) write(*,*) trim(rname)//": callsoil = ",callsoil
440         
441     if (is_master) write(*,*)trim(rname)//&
442       ": Rad transfer is computed every iradia", &
443       " physical timestep"
444     iradia=1 ! default value
445     call getin_p("iradia",iradia)
446     if (is_master) write(*,*)trim(rname)//": iradia = ",iradia
447       
448     if (is_master) write(*,*)trim(rname)//": Rayleigh scattering ?"
449     rayleigh=.false.
450     call getin_p("rayleigh",rayleigh)
451     if (is_master) write(*,*)trim(rname)//": rayleigh = ",rayleigh
452     
453     if (is_master) write(*,*) trim(rname)//&
454       ": prohibit rayleigh calculations outside wavenumber grid?"
455     strictboundrayleigh=.true. ! default value
456     call getin_p("strictboundrayleigh",strictboundrayleigh)
457     if (is_master) write(*,*) trim(rname)//&
458       ": strictboundrayleigh = ",strictboundrayleigh
459
460     if (is_master) write(*,*) trim(rname)//&
461       ": Use blackbody for stellar spectrum ?"
462     stelbbody=.false. ! default value
463     call getin_p("stelbbody",stelbbody)
464     if (is_master) write(*,*) trim(rname)//": stelbbody = ",stelbbody
465
466     if (is_master) write(*,*) trim(rname)//": Stellar blackbody temperature ?"
467     stelTbb=5800.0 ! default value
468     call getin_p("stelTbb",stelTbb)
469     if (is_master) write(*,*) trim(rname)//": stelTbb = ",stelTbb
470
471     if (is_master) write(*,*) trim(rname)//": Output mean OLR in 1D?"
472     meanOLR=.false.
473     call getin_p("meanOLR",meanOLR)
474     if (is_master) write(*,*) trim(rname)//": meanOLR = ",meanOLR
475
476     if (is_master) write(*,*)trim(rname)//": Output spectral OLR in 3D?"
477     specOLR=.false.
478     call getin_p("specOLR",specOLR)
479     if (is_master) write(*,*)trim(rname)//": specOLR = ",specOLR
480
481     if (is_master) write(*,*)trim(rname)//&
482       ": Output diagnostic optical thickness attenuation exp(-tau)?"
483     diagdtau=.false.
484     call getin_p("diagdtau",diagdtau)
485     if (is_master) write(*,*)trim(rname)//": diagdtau = ",diagdtau
486
487     if (is_master) write(*,*)trim(rname)//": Operate in kastprof mode?"
488     kastprof=.false.
489     call getin_p("kastprof",kastprof)
490     if (is_master) write(*,*)trim(rname)//": kastprof = ",kastprof
491
492     if (is_master) write(*,*)trim(rname)//&
493       ": Uniform absorption in radiative transfer?"
494     graybody=.false.
495     call getin_p("graybody",graybody)
496     if (is_master) write(*,*)trim(rname)//": graybody = ",graybody
497
498! Soil model
499     if (is_master) write(*,*)trim(rname)//&
500       ": Number of sub-surface layers for soil scheme?"
501     ! default value of nsoilmx set in comsoil_h
502     call getin_p("nsoilmx",nsoilmx)
503     if (is_master) write(*,*)trim(rname)//": nsoilmx=",nsoilmx
504     
505     if (is_master) write(*,*)trim(rname)//&
506       ": Thickness of topmost soil layer (m)?"
507     ! default value of lay1_soil set in comsoil_h
508     call getin_p("lay1_soil",lay1_soil)
509     if (is_master) write(*,*)trim(rname)//": lay1_soil = ",lay1_soil
510     
511     if (is_master) write(*,*)trim(rname)//&
512       ": Coefficient for soil layer thickness distribution?"
513     ! default value of alpha_soil set in comsoil_h
514     call getin_p("alpha_soil",alpha_soil)
515     if (is_master) write(*,*)trim(rname)//": alpha_soil = ",alpha_soil
516
517! Slab Ocean
518     if (is_master) write(*,*) trim(rname)//": Use slab-ocean ?"
519     ok_slab_ocean=.false.         ! default value
520     call getin_p("ok_slab_ocean",ok_slab_ocean)
521     if (is_master) write(*,*) trim(rname)//": ok_slab_ocean = ",ok_slab_ocean
522     ! Sanity check: for now slab oncean only works in serial mode
523!     if (ok_slab_ocean.and.is_parallel) then
524!       call abort_physic(rname," Error: slab ocean should only be used in serial mode!",1)
525!     endif
526
527!     if (is_master) write(*,*) trim(rname)//": Use slab-sea-ice ?"
528!     ok_slab_sic=.true.         ! default value
529!     call getin_p("ok_slab_sic",ok_slab_sic)
530!     if (is_master) write(*,*) trim(rname)//": ok_slab_sic = ",ok_slab_sic
531
532!     if (is_master) write(*,*) trim(rname)//&
533!       ": Use heat transport for the ocean ?"
534!     ok_slab_heat_transp=.true.   ! default value
535!     call getin_p("ok_slab_heat_transp",ok_slab_heat_transp)
536!     if (is_master) write(*,*) trim(rname)//&
537!       ": ok_slab_heat_transp = ",ok_slab_heat_transp
538
539! Volcanoes (as point sources of tracers)
540      if (is_master) write(*,*) trim(rname)//": Erupt volcano ?"
541      callvolcano=.false. ! default value
542      call getin_p("callvolcano",callvolcano)
543      if (is_master) write(*,*) trim(rname)//": callvolcano = ",callvolcano
544     
545      ! Sanity check: volcano requires tracers
546      if (callvolcano.and.(.not.tracer)) then
547        call abort_physic(rname,&
548             " Error: volcano option requires using tracers",1)
549      endif
550
551! Photochemistry and chemistry in the thermosphere
552
553     if (is_master) write(*,*) trim(rname)//": Use photochemistry ?"
554     photochem=.false.         ! default value
555     call getin_p("photochem",photochem)
556     if (is_master) write(*,*) trim(rname)//": photochem = ",photochem
557
558     if (is_master) write(*,*) trim(rname)//": Use photolysis heat table ?"
559     photoheat=.false.         ! default value
560     call getin_p("photoheat",photoheat)
561     if (is_master) write(*,*) trim(rname)//": photoheat = ",photoheat
562
563     if (is_master) write(*,*) trim(rname)//&
564       ": Use photolysis online calculation ?"
565     jonline=.false.         ! default value
566     call getin_p("jonline",jonline)
567     if (is_master) write(*,*) trim(rname)//": jonline = ",jonline
568
569     if (is_master) write(*,*) trim(rname)//": Use deposition ?"
570     depos=.false.         ! default value
571     call getin_p("depos",depos)
572     if (is_master) write(*,*) trim(rname)//": depos = ",depos
573
574     if (is_master) write(*,*)trim(rname)//": Production of haze ?"
575     haze=.false. ! default value
576     call getin_p("haze",haze)
577     if (is_master) write(*,*)trim(rname)//": haze = ",haze
578
579     if (is_master) write(*,*)trim(rname)//": Output photochemical surface UV flux ?"
580     output_writediagspecUV=.false. ! default value
581     call getin_p("output_writediagspecUV",output_writediagspecUV)
582     if (is_master) write(*,*)trim(rname)//": output_writediagspecUV = ",output_writediagspecUV
583
584     if (is_master) write(*,*)trim(rname)//": Chemical network file ?"
585     reactfile = 'network.def' ! default value
586     call getin_p("reactfile",reactfile)
587     if (is_master) write(*,*)trim(rname)//": chemical network file = ",trim(reactfile)
588
589! Global1D mean and solar zenith angle
590
591     if (ngrid.eq.1) then
592      PRINT*, 'Simulate global averaged conditions ?'
593      global1d = .false. ! default value
594      call getin_p("global1d",global1d)
595      write(*,*) "global1d = ",global1d
596     
597      ! Test of incompatibility : if global1d is true, there should not be any diurnal cycle.
598      if (global1d.and.diurnal) then
599         call abort_physic(rname,'if global1d is true, diurnal must be set to false',1)
600      endif
601
602      if (global1d) then
603         PRINT *,'Solar Zenith angle (deg.) ?'
604         PRINT *,'(assumed for averaged solar flux S/4)'
605         szangle=60.0  ! default value
606         call getin_p("szangle",szangle)
607         write(*,*) "szangle = ",szangle
608      endif
609   endif ! of if (ngrid.eq.1)
610
611! Test of incompatibility:
612! if kastprof used, we must be in 1D
613     if (kastprof.and.(ngrid.gt.1)) then
614       call abort_physic(rname,'kastprof can only be used in 1D!',1)
615     endif
616
617     if (is_master) write(*,*)trim(rname)//&
618       ": Stratospheric temperature for kastprof mode?"
619     Tstrat=167.0
620     call getin_p("Tstrat",Tstrat)
621     if (is_master) write(*,*)trim(rname)//": Tstrat = ",Tstrat
622
623     if (is_master) write(*,*)trim(rname)//": Remove lower boundary?"
624     nosurf=.false.
625     call getin_p("nosurf",nosurf)
626     if (is_master) write(*,*)trim(rname)//": nosurf = ",nosurf
627
628! Tests of incompatibility:
629     if (nosurf.and.callsoil) then
630       if (is_master) then
631         write(*,*)trim(rname)//'nosurf not compatible with soil scheme!'
632         write(*,*)trim(rname)//'... got to make a choice!'
633       endif
634       call abort_physic(rname,"nosurf not compatible with soil scheme!",1)
635     endif
636
637     if (is_master) write(*,*)trim(rname)//": Add an internal heat flux?", &
638                   "... matters only if callsoil=F"
639     intheat=0.
640     call getin_p("intheat",intheat)
641     if (is_master) write(*,*)trim(rname)//": intheat = ",intheat
642
643     if (is_master) write(*,*)trim(rname)//&
644       ": Use Newtonian cooling for radiative transfer?"
645     newtonian=.false.
646     call getin_p("newtonian",newtonian)
647     if (is_master) write(*,*)trim(rname)//": newtonian = ",newtonian
648
649! Tests of incompatibility:
650     if (newtonian.and.corrk) then
651       call abort_physic(rname,'newtonian not compatible with correlated-k!',1)
652     endif
653     if (newtonian.and.calladj) then
654       call abort_physic(rname,'newtonian not compatible with adjustment!',1)
655     endif
656     if (newtonian.and.calldifv) then
657       call abort_physic(rname,'newtonian not compatible with a boundary layer!',1)
658     endif
659
660     if (is_master) write(*,*)trim(rname)//": Test physics timescale in 1D?"
661     testradtimes=.false.
662     call getin_p("testradtimes",testradtimes)
663     if (is_master) write(*,*)trim(rname)//": testradtimes = ",testradtimes
664
665! Test of incompatibility:
666! if testradtimes used, we must be in 1D
667     if (testradtimes.and.(ngrid.gt.1)) then
668       call abort_physic(rname,'testradtimes can only be used in 1D!',1)
669     endif
670
671     if (is_master) write(*,*)trim(rname)//": Default planetary temperature?"
672     tplanet=215.0
673     call getin_p("tplanet",tplanet)
674     if (is_master) write(*,*)trim(rname)//": tplanet = ",tplanet
675
676     if (is_master) write(*,*)trim(rname)//": Which star?"
677     startype=1 ! default value = Sol
678     call getin_p("startype",startype)
679     if (is_master) write(*,*)trim(rname)//": startype = ",startype
680
681     if (is_master) write(*,*)trim(rname)//": Value of stellar flux at 1 AU?"
682     Fat1AU=1356.0 ! default value = Sol today
683     call getin_p("Fat1AU",Fat1AU)
684     if (is_master) write(*,*)trim(rname)//": Fat1AU = ",Fat1AU
685
686
687! TRACERS:
688
689     if (is_master) write(*,*)trim(rname)//": Varying H2O cloud fraction?"
690     CLFvarying=.false.     ! default value
691     call getin_p("CLFvarying",CLFvarying)
692     if (is_master) write(*,*)trim(rname)//": CLFvarying = ",CLFvarying
693
694     if (is_master) write(*,*)trim(rname)//&
695       ": Value of fixed H2O cloud fraction?"
696     CLFfixval=1.0                ! default value
697     call getin_p("CLFfixval",CLFfixval)
698     if (is_master) write(*,*)trim(rname)//": CLFfixval = ",CLFfixval
699
700     if (is_master) write(*,*)trim(rname)//": fixed radii for Cloud particles?"
701     radfixed=.false. ! default value
702     call getin_p("radfixed",radfixed)
703     if (is_master) write(*,*)trim(rname)//": radfixed = ",radfixed
704
705     if(kastprof)then
706        radfixed=.true.
707     endif 
708
709     if (is_master) write(*,*)trim(rname)//&
710       ": Number mixing ratio of CO2 ice particles:"
711     Nmix_co2=1.e6 ! default value
712     call getin_p("Nmix_co2",Nmix_co2)
713     if (is_master) write(*,*)trim(rname)//": Nmix_co2 = ",Nmix_co2
714
715     if (is_master) write(*,*)trim(rname)//&
716       "Number of radiatively active aerosols:"
717     naerkind=0 ! default value
718     call getin_p("naerkind",naerkind)
719     if (is_master) write(*,*)trim(rname)//": naerkind = ",naerkind
720
721     if (is_master) write(*,*)trim(rname)//": Opacity of dust (if used):"
722     dusttau=0. ! default value
723     call getin_p("dusttau",dusttau)
724     if (is_master) write(*,*)trim(rname)//": dusttau = ",dusttau
725
726     if (is_master) write(*,*)trim(rname)//": Radiatively active CO2 aerosols?"
727     aeroco2=.false.     ! default value
728     call getin_p("aeroco2",aeroco2)
729     if (is_master) write(*,*)trim(rname)//": aeroco2 = ",aeroco2
730
731     if (is_master) write(*,*)trim(rname)//": Fixed CO2 aerosol distribution?"
732     aerofixco2=.false.     ! default value
733     call getin_p("aerofixco2",aerofixco2)
734     if (is_master) write(*,*)trim(rname)//": aerofixco2 = ",aerofixco2
735
736     if (is_master) write(*,*)trim(rname)//": Radiatively active water ice?"
737     aeroh2o=.false.     ! default value
738     call getin_p("aeroh2o",aeroh2o)
739     if (is_master) write(*,*)trim(rname)//": aeroh2o = ",aeroh2o
740
741     if (is_master) write(*,*)trim(rname)//": Fixed H2O aerosol distribution?"
742     aerofixh2o=.false.     ! default value
743     call getin_p("aerofixh2o",aerofixh2o)
744     if (is_master) write(*,*)trim(rname)//": aerofixh2o = ",aerofixh2o
745
746     if (is_master) write(*,*)trim(rname)//&
747       ": Radiatively active sulfuric acid aerosols?"
748     aeroh2so4=.false.     ! default value
749     call getin_p("aeroh2so4",aeroh2so4)
750     if (is_master) write(*,*)trim(rname)//": aeroh2so4 = ",aeroh2so4
751 
752     if (is_master) write(*,*)trim(rname)//&
753       ": Radiatively active auroral aerosols?"
754     aeroaurora=.false.     ! default value
755     call getin_p("aeroaurora",aeroaurora)
756     if (is_master) write(*,*)trim(rname)//": aeroaurora = ",aeroaurora
757
758     if (is_master) write(*,*)trim(rname)//&
759     ": Radiatively active Generic Condensable aerosols?"
760     aerogeneric=0     ! default value
761     call getin_p("aerogeneric",aerogeneric)
762     if (is_master) write(*,*)trim(rname)//":aerogeneric",aerogeneric
763
764
765!=================================
766! TWOLAY scheme and NH3 cloudare left for retrocompatibility only,
767! You should now use N-LAYER scheme (see below).
768
769     if (is_master) write(*,*)trim(rname)//&
770       ": Radiatively active two-layer aerosols?"
771     aeroback2lay=.false.     ! default value
772     call getin_p("aeroback2lay",aeroback2lay)
773     if (is_master) write(*,*)trim(rname)//": aeroback2lay = ",aeroback2lay
774
775     if (aeroback2lay.and.is_master) then
776       print*,'Warning : The TWOLAY AEROSOL scheme is deprecated and buggy...'
777       print*,'You should use the generic n-layer scheme (see aeronlay).'
778     endif
779
780     if (is_master) write(*,*)trim(rname)//": Radiatively active ammonia cloud?"
781     aeronh3=.false.     ! default value
782     call getin_p("aeronh3",aeronh3)
783     if (is_master) write(*,*)trim(rname)//": aeronh3 = ",aeronh3
784
785     if (aeronh3.and.is_master) then
786       print*,'Warning : You are using specific NH3 cloud scheme ...'
787       print*,'You should use the generic n-layer scheme (see aeronlay).'
788     endif
789
790     if (is_master) write(*,*)"Radiatively active Venus clouds ?"
791     aerovenus=.false. ! default value
792     call getin_p("aerovenus",aerovenus)
793     if (aerovenus) then
794       aerovenus1=.true.     ! default value
795       aerovenus2=.true.     ! default value
796       aerovenus2p=.true.     ! default value
797       aerovenus3=.true.     ! default value
798       aerovenusUV=.true.     ! default value
799     else
800       aerovenus1=.false.     ! default value
801       aerovenus2=.false.     ! default value
802       aerovenus2p=.false.     ! default value
803       aerovenus3=.false.     ! default value
804       aerovenusUV=.false.     ! default value
805     endif
806     ! in case the user wants to specifically set/unset sub-options
807     call getin_p("aerovenus1",aerovenus1)
808     if (is_master) write(*,*)" aerovenus1 = ",aerovenus1
809     call getin_p("aerovenus2",aerovenus2)
810     if (is_master) write(*,*)" aerovenus2 = ",aerovenus2
811     call getin_p("aerovenus2p",aerovenus2p)
812     if (is_master) write(*,*)" aerovenus2p= ",aerovenus2p
813     call getin_p("aerovenus3",aerovenus3)
814     if (is_master) write(*,*)" aerovenus3 = ",aerovenus3
815     call getin_p("aerovenusUV",aerovenusUV)
816     if (is_master) write(*,*)" aerovenusUV= ",aerovenusUV
817     ! Sanity check: if any of the aerovenus* is set to true
818     ! then aeronovenus should also be set to true
819     if ((.not.aerovenus).and.(aerovenus1.or.aerovenus2.or.aerovenus2p.or.&
820                               aerovenus3.or.aerovenusUV)) then
821      if(is_master) then
822       write(*,*)" Error, if you set some of the aerovenus* to true"
823       write(*,*)" then flag aerovenus should be set to true as well!"
824      endif
825      call abort_physic("inifis"," aerovenus* flags mismatch!",1)
826     endif
827     
828     if (is_master) write(*,*)trim(rname)//&
829       ": TWOLAY AEROSOL: total optical depth "//&
830       "in the tropospheric layer (visible)"
831     obs_tau_col_tropo=8.D0
832     call getin_p("obs_tau_col_tropo",obs_tau_col_tropo)
833     if (is_master) write(*,*)trim(rname)//&
834       ": obs_tau_col_tropo = ",obs_tau_col_tropo
835
836     if (is_master) write(*,*)trim(rname)//&
837       ": TWOLAY AEROSOL: total optical depth "//&
838       "in the stratospheric layer (visible)"
839     obs_tau_col_strato=0.08D0
840     call getin_p("obs_tau_col_strato",obs_tau_col_strato)
841     if (is_master) write(*,*)trim(rname)//&
842       ": obs_tau_col_strato = ",obs_tau_col_strato
843
844     if (is_master) write(*,*)trim(rname)//&
845       ": TWOLAY AEROSOL: optprop_back2lay_vis?"
846     optprop_back2lay_vis = 'optprop_saturn_vis_n20.dat'
847     call getin_p("optprop_back2lay_vis",optprop_back2lay_vis)
848     if (is_master) write(*,*)trim(rname)//&
849       ": optprop_back2lay_vis = ",trim(optprop_back2lay_vis)
850
851     if (is_master) write(*,*)trim(rname)//&
852       ": TWOLAY AEROSOL: optprop_back2lay_ir?"
853     optprop_back2lay_ir = 'optprop_saturn_ir_n20.dat'
854     call getin_p("optprop_back2lay_ir",optprop_back2lay_ir)
855     if (is_master) write(*,*)trim(rname)//&
856       ": optprop_back2lay_ir = ",trim(optprop_back2lay_ir)
857     
858     if (is_master) write(*,*)trim(rname)//&
859       ": TWOLAY AEROSOL: pres_bottom_tropo? in pa"
860     pres_bottom_tropo=66000.0
861     call getin_p("pres_bottom_tropo",pres_bottom_tropo)
862     if (is_master) write(*,*)trim(rname)//&
863       ": pres_bottom_tropo = ",pres_bottom_tropo
864
865     if (is_master) write(*,*)trim(rname)//&
866       ": TWOLAY AEROSOL: pres_top_tropo? in pa"
867     pres_top_tropo=18000.0
868     call getin_p("pres_top_tropo",pres_top_tropo)
869     if (is_master) write(*,*)trim(rname)//&
870       ": pres_top_tropo = ",pres_top_tropo
871
872     if (is_master) write(*,*)trim(rname)//&
873       ": TWOLAY AEROSOL: pres_bottom_strato? in pa"
874     pres_bottom_strato=2000.0
875     call getin_p("pres_bottom_strato",pres_bottom_strato)
876     if (is_master) write(*,*)trim(rname)//&
877       ": pres_bottom_strato = ",pres_bottom_strato
878
879     ! Sanity check
880     if (pres_bottom_strato .gt. pres_top_tropo) then
881       if(is_master) then
882       print*,'Error : TWOLAY AEROSOL, Please ensure that in callphys.def'
883       print*,'you have pres_top_tropo > pres_bottom_strato !'
884       endif
885       call abort_physic(rname," pres_top_tropo > pres_bottom_strato!",1)
886     endif
887
888     if (is_master) write(*,*)trim(rname)//&
889       ": TWOLAY AEROSOL: pres_top_strato? in pa"
890     pres_top_strato=100.0
891     call getin_p("pres_top_strato",pres_top_strato)
892     if (is_master) write(*,*)trim(rname)//&
893       ": pres_top_strato = ",pres_top_strato
894
895     if (is_master) write(*,*)trim(rname)//&
896       ": TWOLAY AEROSOL: particle size in the ", &
897       "tropospheric layer, in meters"
898     size_tropo=2.e-6
899     call getin_p("size_tropo",size_tropo)
900     if (is_master) write(*,*)trim(rname)//": size_tropo = ",size_tropo
901
902     if (is_master) write(*,*)trim(rname)//&
903       ": TWOLAY AEROSOL: particle size in the ", &
904       "stratospheric layer, in meters"
905     size_strato=1.e-7
906     call getin_p("size_strato",size_strato)
907     if (is_master) write(*,*)trim(rname)//": size_strato = ",size_strato
908
909     if (is_master) write(*,*)trim(rname)//&
910       ": NH3 (thin) cloud: total optical depth"
911     tau_nh3_cloud=7.D0
912     call getin_p("tau_nh3_cloud",tau_nh3_cloud)
913     if (is_master) write(*,*)trim(rname)//": tau_nh3_cloud = ",tau_nh3_cloud
914
915     if (is_master) write(*,*)trim(rname)//": NH3 (thin) cloud pressure level"
916     pres_nh3_cloud=7.D0
917     call getin_p("pres_nh3_cloud",pres_nh3_cloud)
918     if (is_master) write(*,*)trim(rname)//": pres_nh3_cloud = ",pres_nh3_cloud
919
920     if (is_master) write(*,*)trim(rname)//": NH3 (thin) cloud: particle sizes"
921     size_nh3_cloud=3.e-6
922     call getin_p("size_nh3_cloud",size_nh3_cloud)
923     if (is_master) write(*,*)trim(rname)//": size_nh3_cloud = ",size_nh3_cloud
924
925!=================================
926! Generic N-LAYER aerosol scheme
927
928     if (is_master) write(*,*)trim(rname)//&
929       ": Radiatively active generic n-layer aerosols?"
930     aeronlay=.false.     ! default value
931     call getin_p("aeronlay",aeronlay)
932     if (is_master) write(*,*)trim(rname)//": aeronlay = ",aeronlay
933
934     if (is_master) write(*,*)trim(rname)//&
935       ": Number of generic aerosols layers?"
936     nlayaero=1     ! default value
937     call getin_p("nlayaero",nlayaero)
938     ! Avoid to allocate arrays of size 0
939     if (aeronlay .and. nlayaero.lt.1) then
940       if (is_master) then
941       print*, " You are trying to set no generic aerosols..."
942       print*, " Set aeronlay=.false. instead ! I abort."
943       endif
944       call abort_physic(rname,"no generic aerosols...",1)
945     endif
946     if (.not. aeronlay) nlayaero=1
947     if (is_master) write(*,*)trim(rname)//": nlayaero = ",nlayaero
948
949     ! This is necessary, we just set the number of aerosol layers
950     IF(.NOT.ALLOCATED(aeronlay_tauref))      ALLOCATE(aeronlay_tauref(nlayaero))     
951     IF(.NOT.ALLOCATED(aeronlay_lamref))      ALLOCATE(aeronlay_lamref(nlayaero))     
952     IF(.NOT.ALLOCATED(aeronlay_choice))      ALLOCATE(aeronlay_choice(nlayaero))     
953     IF(.NOT.ALLOCATED(aeronlay_pbot))        ALLOCATE(aeronlay_pbot(nlayaero))     
954     IF(.NOT.ALLOCATED(aeronlay_ptop))        ALLOCATE(aeronlay_ptop(nlayaero))     
955     IF(.NOT.ALLOCATED(aeronlay_sclhght))     ALLOCATE(aeronlay_sclhght(nlayaero))     
956     IF(.NOT.ALLOCATED(aeronlay_size))        ALLOCATE(aeronlay_size(nlayaero))     
957     IF(.NOT.ALLOCATED(aeronlay_nueff))       ALLOCATE(aeronlay_nueff(nlayaero))     
958     IF(.NOT.ALLOCATED(optprop_aeronlay_ir))  ALLOCATE(optprop_aeronlay_ir(nlayaero))     
959     IF(.NOT.ALLOCATED(optprop_aeronlay_vis)) ALLOCATE(optprop_aeronlay_vis(nlayaero))     
960
961     if (is_master) write(*,*)trim(rname)//&
962       ": Generic n-layer aerosols: Optical depth at reference wavelenght"
963     aeronlay_tauref=1.0E-1
964     call getin_p("aeronlay_tauref",aeronlay_tauref)
965     if (is_master) write(*,*)trim(rname)//&
966       ": aeronlay_tauref = ",aeronlay_tauref
967
968     if (is_master) write(*,*)trim(rname)//&
969       ": Generic n-layer aerosols: Reference wavelenght for optical depths (m)"
970     aeronlay_lamref=0.6E-6
971     call getin_p("aeronlay_lamref",aeronlay_lamref)
972     if (is_master) write(*,*)trim(rname)//&
973       ": aeronlay_lamref = ",aeronlay_lamref
974
975     if (is_master) then
976       write(*,*)trim(rname)//&
977       ": Generic n-layer aerosols: Vertical profile choice : "
978       write(*,*)trim(rname)//&
979       "             (1) Tau btwn ptop and pbot follows atm. scale height"
980       write(*,*)trim(rname)//&
981       "         or  (2) Tau above pbot follows its own scale height"
982     endif
983     aeronlay_choice=1
984     call getin_p("aeronlay_choice",aeronlay_choice)
985     if (is_master) write(*,*)trim(rname)//&
986       ": aeronlay_choice = ",aeronlay_choice
987
988     if (is_master) write(*,*)trim(rname)//&
989       ": Generic n-layer aerosols: bottom pressures (Pa)"
990     aeronlay_pbot=2000.0
991     call getin_p("aeronlay_pbot",aeronlay_pbot)
992     if (is_master) write(*,*)trim(rname)//": aeronlay_pbot = ",aeronlay_pbot
993
994     if (is_master) write(*,*)trim(rname)//&
995       ": Generic n-layer aerosols: (if choice=1) Top pressures (Pa) "
996     aeronlay_ptop=300000.0
997     call getin_p("aeronlay_ptop",aeronlay_ptop)
998     if (is_master) write(*,*)trim(rname)//": aeronlay_ptop = ",aeronlay_ptop
999
1000     if (is_master) write(*,*)trim(rname)//&
1001       ": Generic n-layer aerosols: (if choice=2) Scale height / atm. scale height"
1002     aeronlay_sclhght=0.2
1003     call getin_p("aeronlay_sclhght",aeronlay_sclhght)
1004     if (is_master) write(*,*)trim(rname)//&
1005       ": aeronlay_sclhght = ",aeronlay_sclhght
1006
1007     if (is_master) write(*,*)trim(rname)//&
1008       ": Generic n-layer aerosols: particles effective radii (m)"
1009     aeronlay_size=1.e-6
1010     call getin_p("aeronlay_size",aeronlay_size)
1011     if (is_master) write(*,*)trim(rname)//": aeronlay_size = ",aeronlay_size
1012
1013     if (is_master) write(*,*)trim(rname)//&
1014       ": Generic n-layer aerosols: particles radii effective variance"
1015     aeronlay_nueff=0.1
1016     call getin_p("aeronlay_nueff",aeronlay_nueff)
1017     if (is_master) write(*,*)trim(rname)//": aeronlay_nueff = ",aeronlay_nueff
1018
1019     if (is_master) write(*,*)trim(rname)//&
1020       ": Generic n-layer aerosols: VIS optical properties file"
1021     optprop_aeronlay_vis = 'optprop_saturn_vis_n20.dat'
1022     call getin_p("optprop_aeronlay_vis",optprop_aeronlay_vis)
1023     if (is_master) write(*,*)trim(rname)//&
1024       ": optprop_aeronlay_vis = ",optprop_aeronlay_vis
1025
1026     if (is_master) write(*,*)trim(rname)//&
1027     ": Generic n-layer aerosols: IR optical properties file"
1028     optprop_aeronlay_ir = 'optprop_saturn_ir_n20.dat'
1029     call getin_p("optprop_aeronlay_ir",optprop_aeronlay_ir)
1030     if (is_master) write(*,*)trim(rname)//&
1031       ": optprop_aeronlay_ir = ",optprop_aeronlay_ir
1032     
1033
1034!=================================
1035
1036     if (is_master) write(*,*)trim(rname)//&
1037       ": Cloud pressure level (with kastprof only):"
1038     cloudlvl=0. ! default value
1039     call getin_p("cloudlvl",cloudlvl)
1040     if (is_master) write(*,*)trim(rname)//": cloudlvl = ",cloudlvl
1041
1042     if (is_master) write(*,*)trim(rname)//&
1043       ": Is the variable gas species radiatively active?"
1044     Tstrat=167.0
1045     varactive=.false.
1046     call getin_p("varactive",varactive)
1047     if (is_master) write(*,*)trim(rname)//": varactive = ",varactive
1048
1049     if (is_master) write(*,*)trim(rname)//&
1050       ": Is the variable gas species distribution set?"
1051     varfixed=.false.
1052     call getin_p("varfixed",varfixed)
1053     if (is_master) write(*,*)trim(rname)//": varfixed = ",varfixed
1054
1055     if (is_master) write(*,*)trim(rname)//&
1056       ": What is the saturation % of the variable species?"
1057     satval=0.8
1058     call getin_p("satval",satval)
1059     if (is_master) write(*,*)trim(rname)//": satval = ",satval
1060
1061
1062! Test of incompatibility:
1063! if varactive, then varfixed should be false
1064     if (varactive.and.varfixed) then
1065       call abort_physic(rname,'if varactive, varfixed must be OFF!',1)
1066     endif
1067
1068     if (is_master) write(*,*)trim(rname)//": Gravitationnal sedimentation ?"
1069     sedimentation=.false. ! default value
1070     call getin_p("sedimentation",sedimentation)
1071     if (is_master) write(*,*)trim(rname)//": sedimentation = ",sedimentation
1072
1073     if (is_master) write(*,*)trim(rname)//": Generic Condensation of tracers ?"
1074     generic_condensation=.false. !default value
1075     call getin_p("generic_condensation",generic_condensation)
1076     if (is_master) write(*,*)trim(rname)//": generic_condensation = ",generic_condensation
1077     
1078     if (is_master) write(*,*)trim(rname)//": Generic rain of tracers ?"
1079     generic_rain=.false. !default value
1080     call getin_p("generic_rain",generic_rain)
1081     if (is_master) write(*,*)trim(rname)//": generic_rain = ",generic_rain
1082
1083     if (is_master) write(*,*)trim(rname)//": Include moist adjustement for GCS ?"
1084     moistadjustment_generic=.false. ! default value
1085     call getin_p("moistadjustment_generic",moistadjustment_generic)
1086     if (is_master) write(*,*)trim(rname)//": moistadjustment_generic = ", moistadjustment_generic
1087
1088     if (is_master) write(*,*)trim(rname)//": Moist convection inhibition for GCS ?"
1089     moist_convection_inhibition=.false. ! default value
1090     call getin_p("moist_convection_inhibition",moist_convection_inhibition)
1091     if (is_master) write(*,*)trim(rname)//": moist_convection_inhibition = ", moist_convection_inhibition
1092     
1093     if (is_master) write(*,*)trim(rname)//": Virtual theta correction ?"
1094     virtual_theta_correction=.false. !default value
1095     call getin_p("virtual_theta_correction",virtual_theta_correction)
1096     if (is_master) write(*,*)trim(rname)//": virtual_theta_correction = ",virtual_theta_correction
1097
1098     if (is_master) write(*,*)trim(rname)//": Compute water cycle ?"
1099     water=.false. ! default value
1100     call getin_p("water",water)
1101     if (is_master) write(*,*)trim(rname)//": water = ",water
1102         
1103! Test of incompatibility:
1104! if water is true, there should be at least a tracer
1105     if (water.and.(.not.tracer)) then
1106        call abort_physic(rname,'if water is ON, tracer must be ON too!',1)
1107     endif
1108
1109     if (is_master) write(*,*)trim(rname)//": Include water condensation ?"
1110     watercond=.false. ! default value
1111     call getin_p("watercond",watercond)
1112     if (is_master) write(*,*)trim(rname)//": watercond = ",watercond
1113
1114! Test of incompatibility:
1115! if watercond is used, then water should be used too
1116     if (watercond.and.(.not.water)) then
1117        call abort_physic(rname,&
1118          'if watercond is used, water should be used too',1)
1119     endif
1120
1121     if (is_master) write(*,*)trim(rname)//": Include water moist adjustement ?"
1122     moistadjustment=.true. ! default value
1123     call getin_p("moistadjustment",moistadjustment)
1124     if (is_master) write(*,*)trim(rname)//": moistadjustment = ", moistadjustment
1125
1126     if (is_master) write(*,*)trim(rname)//": Include water precipitation ?"
1127     waterrain=.false. ! default value
1128     call getin_p("waterrain",waterrain)
1129     if (is_master) write(*,*)trim(rname)//": waterrain = ",waterrain
1130
1131     if (is_master) write(*,*)trim(rname)//": Include surface hydrology ?"
1132     hydrology=.false. ! default value
1133     call getin_p("hydrology",hydrology)
1134     if (is_master) write(*,*)trim(rname)//": hydrology = ",hydrology
1135
1136     if (is_master) write(*,*)trim(rname)//": Spectral Dependant albedo ?"
1137     albedo_spectral_mode=.false. ! default value
1138     call getin_p("albedo_spectral_mode",albedo_spectral_mode)
1139     if (is_master) write(*,*)trim(rname)//&
1140     ": albedo_spectral_mode = ",albedo_spectral_mode
1141
1142     if (is_master) then
1143       write(*,*)trim(rname)//": Snow albedo ?"
1144       write(*,*)trim(rname)//": If albedo_spectral_mode=.true., then this "
1145       write(*,*)trim(rname)//": corresponds to the 0.5 microns snow albedo."
1146     endif
1147     albedosnow=0.5         ! default value
1148     call getin_p("albedosnow",albedosnow)
1149     if (is_master) write(*,*)trim(rname)//": albedosnow = ",albedosnow
1150         
1151     if (is_master) write(*,*)trim(rname)//": Ocean albedo ?"
1152     alb_ocean=0.07         ! default value
1153     call getin_p("alb_ocean",alb_ocean)
1154     if (is_master) write(*,*)trim(rname)//": alb_ocean = ",alb_ocean
1155         
1156     if (is_master) write(*,*)trim(rname)//": CO2 ice albedo ?"
1157     albedoco2ice=0.5       ! default value
1158     call getin_p("albedoco2ice",albedoco2ice)
1159     if (is_master) write(*,*)trim(rname)//": albedoco2ice = ",albedoco2ice
1160
1161     if (is_master) write(*,*)trim(rname)//": Maximum ice thickness ?"
1162     maxicethick=2.0         ! default value
1163     call getin_p("maxicethick",maxicethick)
1164     if (is_master) write(*,*)trim(rname)//": maxicethick = ",maxicethick
1165
1166     if (is_master) write(*,*)trim(rname)//": Freezing point of seawater ?"
1167     Tsaldiff=-1.8          ! default value
1168     call getin_p("Tsaldiff",Tsaldiff)
1169     if (is_master) write(*,*)trim(rname)//": Tsaldiff = ",Tsaldiff
1170
1171     if (is_master) write(*,*)trim(rname)//": Minimum eddy mix coeff in 1D ?"
1172     kmixmin=1.0e-2         ! default value
1173     call getin_p("kmixmin",kmixmin)
1174     if (is_master) write(*,*)trim(rname)//": kmixmin = ",kmixmin
1175     
1176     varspec=.false. !default value
1177     call getin_p("varspec",varspec)
1178     if (varspec) then
1179       call getin_p("varspec_data",varspec_data)
1180       call getin_p("nvarlayer",nvarlayer)
1181     endif
1182
1183     if (is_master) write(*,*)'Predefined Cp from dynamics is ',cpp,'J kg^-1 K^-1'
1184     if (is_master) write(*,*)'Predefined Mg from dynamics is ',mugaz,'amu'
1185
1186     force_cpp=.false. ! default value
1187     call getin_p("force_cpp",force_cpp)
1188     if (force_cpp) then
1189      if (is_master) write(*,*)trim(rname)//": force_cpp = ",force_cpp
1190      if (is_master) write(*,*)trim(rname)//": force_cpp is deprecated.",&
1191      "Set cpp_mugaz_mode=1 in callfis to emulate force_cpp=.true."
1192      call abort_physic(rname,"Anyway, you need to set force_cpp=.false. to continue.",1)
1193     endif
1194
1195     if (is_master) write(*,*)trim(rname)//&
1196     ": where do you want your cpp/mugaz value to come from?",&
1197     "=> 0: dynamics (3d), 1: forced in callfis (1d), 2: computed from gases.def (1d)?"
1198     cpp_mugaz_mode = 0 ! default value
1199     call getin_p("cpp_mugaz_mode",cpp_mugaz_mode)
1200     if (is_master) write(*,*)trim(rname)//": cpp_mugaz_mode = ",cpp_mugaz_mode
1201
1202     if ((cpp_mugaz_mode >= 1).and.(is_master).and.(ngrid>1)) then
1203        write(*,*)'    !!!!! Be aware that having different values for cpp and mugaz in the dynamics and physics'
1204        write(*,*)'    in 3D can result in very pathological behavior. You have been warned !!!!!'
1205      else if ((cpp_mugaz_mode == 0).and.(is_master).and.(ngrid==1)) then
1206        ! for this specific 1D error we will remove run.def before aborting JL22
1207        call system("rm -rf run.def")
1208        call abort_physic(rname,"cpp_mugaz_mode must be >= 1 in 1d",1)
1209     endif
1210
1211     if (cpp_mugaz_mode == 1) then
1212       mugaz = -99999.
1213       if (is_master) write(*,*)trim(rname)//&
1214         ": MEAN MOLECULAR MASS in g mol-1 ?"
1215       call getin_p("mugaz",mugaz)
1216       IF (mugaz.eq.-99999.) THEN
1217         call abort_physic(rname,"mugaz must be set if cpp_mugaz_mode = 1",1)
1218       ENDIF
1219       cpp = -99999.
1220       if (is_master) write(*,*)trim(rname)//&
1221         ": SPECIFIC HEAT CAPACITY in J K-1 kg-1 ?"
1222       call getin_p("cpp",cpp)
1223       IF (cpp.eq.-99999.) THEN
1224           PRINT *, "cpp must be set if cpp_mugaz_mode = 1"
1225           STOP
1226       ENDIF
1227       if (is_master) write(*,*)'New Cp from callfis is ',cpp,'J kg^-1 K^-1'
1228       if (is_master) write(*,*)'New Mg from callfis is ',mugaz,'amu'
1229 
1230     endif ! of if (cpp_mugaz_mode == 1)
1231     call su_gases
1232     call calc_cpp_mugaz
1233
1234     if (is_master) then
1235       PRINT*,'--------------------------------------------'
1236       PRINT*
1237       PRINT*
1238     endif
1239  ELSE
1240     call abort_physic(rname,'Cannot read file callphys.def. Is it here ?',1)
1241  ENDIF ! of IF(iscallphys)
1242
1243  if (is_master) then
1244    PRINT*
1245    PRINT*,'inifis: daysec',daysec
1246    PRINT*
1247    PRINT*,'inifis: The radiative transfer is computed:'
1248    PRINT*,'           each ',iradia,' physical time-step'
1249    PRINT*,'        or each ',iradia*dtphys,' seconds'
1250    PRINT*
1251  endif
1252
1253!-----------------------------------------------------------------------
1254!     Some more initialization:
1255!     ------------------------
1256
1257  ! Initializations for comgeomfi_h
1258#ifndef MESOSCALE
1259  totarea=SSUM(ngrid,parea,1)
1260  call planetwide_sumval(parea,totarea_planet)
1261
1262  !! those are defined in comdiurn_h.F90
1263  IF (.not.ALLOCATED(sinlat)) ALLOCATE(sinlat(ngrid))
1264  IF (.not.ALLOCATED(coslat)) ALLOCATE(coslat(ngrid))
1265  IF (.not.ALLOCATED(sinlon)) ALLOCATE(sinlon(ngrid))
1266  IF (.not.ALLOCATED(coslon)) ALLOCATE(coslon(ngrid))
1267
1268  DO ig=1,ngrid
1269     sinlat(ig)=sin(plat(ig))
1270     coslat(ig)=cos(plat(ig))
1271     sinlon(ig)=sin(plon(ig))
1272     coslon(ig)=cos(plon(ig))
1273  ENDDO
1274#endif
1275  ! initialize variables in radinc_h
1276  call ini_radinc_h(nlayer,tplanckmin,tplanckmax,dtplanck)
1277
1278  ! initialize variables and allocate arrays in radcommon_h
1279  call ini_radcommon_h(naerkind)
1280   
1281  ! allocate "comsoil_h" arrays
1282  call ini_comsoil_h(ngrid)
1283   
1284  END SUBROUTINE inifis
1285
1286END MODULE inifis_mod
Note: See TracBrowser for help on using the repository browser.