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

Last change on this file since 1351 was 1315, checked in by milmd, 10 years ago

LMDZ.GENERIC. OpenMP directives added in generic physic. When running in pure OpenMP or hybrid OpenMP/MPI, may have some bugs with condense_cloud and wstats routines.

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