source: trunk/LMDZ.GENERIC/libf/phystd/inifis.F @ 1477

Last change on this file since 1477 was 1423, checked in by sglmd, 10 years ago

Second part of solving the top-of-atmosphere problem in the radiative transfer: Return to a previous interpolation scheme

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