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

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

Bug correction in the initialization of specific gas constant in Saturn physics.
EM

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