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

Last change on this file since 2635 was 2635, checked in by aslmd, 3 years ago

Introduction of cpp_mugaz_mode to decide what source should be used for cpp and mugaz. Force_cpp and check_cpp_match are deprecated (see issue la-communaut-des-mod-les-atmosph-riques-plan-taires/git-trunk#27).

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