source: trunk/LMDZ.MARS/libf/phymars/inifis.F @ 553

Last change on this file since 553 was 553, checked in by acolaitis, 13 years ago

Sorry, forgot a case in inifis.

File size: 26.9 KB
Line 
1      SUBROUTINE inifis(
2     $           ngrid,nlayer
3     $           ,day_ini,pdaysec,ptimestep
4     $           ,plat,plon,parea
5     $           ,prad,pg,pr,pcpp
6#ifdef MESOSCALE
7#include "meso_inc/meso_inc_inifisinvar.F"
8#endif
9     $           )
10!
11!=======================================================================
12!
13!   purpose:
14!   -------
15!
16!   Initialisation for the physical parametrisations of the LMD
17!   martian atmospheric general circulation modele.
18!
19!   author: Frederic Hourdin 15 / 10 /93
20!   -------
21!   modified: Sebastien Lebonnois 11/06/2003 (new callphys.def)
22!             Ehouarn Millour (oct. 2008) tracers are now identified
23!              by their names and may not be contiguously
24!              stored in the q(:,:,:,:) array
25!             E.M. (june 2009) use getin routine to load parameters
26!             adapted to the mesoscale use - Aymeric Spiga - 01/2007-07/2011
27!
28!
29!   arguments:
30!   ----------
31!
32!   input:
33!   ------
34!
35!    ngrid                 Size of the horizontal grid.
36!                          All internal loops are performed on that grid.
37!    nlayer                Number of vertical layers.
38!    pdayref               Day of reference for the simulation
39!    pday                  Number of days counted from the North. Spring
40!                          equinoxe.
41!
42!=======================================================================
43!
44!-----------------------------------------------------------------------
45!   declarations:
46!   -------------
47! to use  'getin'
48      USE ioipsl_getincom
49      IMPLICIT NONE
50#include "dimensions.h"
51#include "dimphys.h"
52#include "planete.h"
53#include "comcstfi.h"
54#include "comsaison.h"
55#include "comdiurn.h"
56#include "comgeomfi.h"
57#include "callkeys.h"
58#include "surfdat.h"
59#include "dimradmars.h"
60#include "yomaer.h"
61#include "datafile.h"
62#include "slope.h"
63#include "microphys.h"
64#include "tracer.h"
65#ifdef MESOSCALE
66#include "comsoil.h"     !!MESOSCALE -- needed to fill volcapa
67#include "meso_inc/meso_inc_inifisvar.F"
68#endif
69      REAL prad,pg,pr,pcpp,pdaysec
70
71      REAL ptimestep
72      INTEGER day_ini
73
74      INTEGER ngrid,nlayer
75      REAL plat(ngrid),plon(ngrid),parea(ngridmx)
76      INTEGER ig,ierr
77 
78!      EXTERNAL iniorbit,orbite
79      EXTERNAL SSUM
80      REAL SSUM
81 
82      CHARACTER ch1*12
83      CHARACTER ch80*80
84
85!      logical chem, h2o
86
87!      chem = .false.
88!      h2o = .false.
89
90      rad=prad
91      cpp=pcpp
92      g=pg
93      r=pr
94      rcp=r/cpp
95      daysec=pdaysec
96      dtphys=ptimestep
97#ifdef MESOSCALE
98#include "meso_inc/meso_inc_inifisini.F"
99#endif
100
101! --------------------------------------------------------
102!     The usual Tests
103!     --------------
104      IF (nlayer.NE.nlayermx) THEN
105         PRINT*,'STOP in inifis'
106         PRINT*,'Probleme de dimensions :'
107         PRINT*,'nlayer     = ',nlayer
108         PRINT*,'nlayermx   = ',nlayermx
109         STOP
110      ENDIF
111
112      IF (ngrid.NE.ngridmx) THEN
113         PRINT*,'STOP in inifis'
114         PRINT*,'Probleme de dimensions :'
115         PRINT*,'ngrid     = ',ngrid
116         PRINT*,'ngridmx   = ',ngridmx
117         STOP
118      ENDIF
119
120! --------------------------------------------------------------
121!  Reading the "callphys.def" file controlling some key options
122! --------------------------------------------------------------
123     
124      ! check that 'callphys.def' file is around
125      OPEN(99,file='callphys.def',status='old',form='formatted'
126     &     ,iostat=ierr)
127      CLOSE(99)
128     
129      IF(ierr.EQ.0) THEN
130         PRINT*
131         PRINT*
132         PRINT*,'--------------------------------------------'
133         PRINT*,' inifis: Parameters for the physics (callphys.def)'
134         PRINT*,'--------------------------------------------'
135
136         write(*,*) "Directory where external input files are:"
137         datafile="/u/forget/WWW/datagcm/datafile"
138         call getin("datadir",datafile) ! default path
139         write(*,*) " datafile = ",trim(datafile)
140
141         write(*,*) "Run with or without tracer transport ?"
142         tracer=.false. ! default value
143         call getin("tracer",tracer)
144         write(*,*) " tracer = ",tracer
145
146         write(*,*) "Diurnal cycle ?"
147         write(*,*) "(if diurnal=False, diurnal averaged solar heating)"
148         diurnal=.true. ! default value
149         call getin("diurnal",diurnal)
150         write(*,*) " diurnal = ",diurnal
151
152         write(*,*) "Seasonal cycle ?"
153         write(*,*) "(if season=False, Ls stays constant, to value ",
154     &   "set in 'start'"
155         season=.true. ! default value
156         call getin("season",season)
157         write(*,*) " season = ",season
158
159         write(*,*) "Write some extra output to the screen ?"
160         lwrite=.false. ! default value
161         call getin("lwrite",lwrite)
162         write(*,*) " lwrite = ",lwrite
163
164         write(*,*) "Save statistics in file stats.nc ?"
165#ifdef MESOSCALE
166         callstats=.false. ! default value
167#else
168         callstats=.true. ! default value
169#endif
170         call getin("callstats",callstats)
171         write(*,*) " callstats = ",callstats
172
173         write(*,*) "Save EOF profiles in file 'profiles' for ",
174     &              "Climate Database?"
175         calleofdump=.false. ! default value
176         call getin("calleofdump",calleofdump)
177         write(*,*) " calleofdump = ",calleofdump
178
179         write(*,*) "Dust scenario: 1=constant dust (read from startfi",
180     &   " or set as tauvis); 2=Viking scenario; =3 MGS scenario,",
181     &   "=4 Mars Year 24 from TES assimilation, ",
182     &   "=24,25 or 26 :Mars Year 24,25 or 26 from TES assimilation"
183         iaervar=3 ! default value
184         call getin("iaervar",iaervar)
185         write(*,*) " iaervar = ",iaervar
186
187         write(*,*) "Reference (visible) dust opacity at 700 Pa ",
188     &   "(matters only if iaervar=1)"
189         ! NB: default value of tauvis is set/read in startfi.nc file
190         call getin("tauvis",tauvis)
191         write(*,*) " tauvis = ",tauvis
192
193         write(*,*) "Dust vertical distribution:"
194         write(*,*) "(=1 top set by topdustref parameter;",
195     & " =2 Viking scenario; =3 MGS scenario)"
196         iddist=3 ! default value
197         call getin("iddist",iddist)
198         write(*,*) " iddist = ",iddist
199
200         write(*,*) "Dust top altitude (km). (Matters only if iddist=1)"
201         topdustref= 90.0 ! default value
202         call getin("topdustref",topdustref)
203         write(*,*) " topdustref = ",topdustref
204
205         write(*,*) "Prescribed surface thermal flux (H/(rho*cp),K m/s)"
206         tke_heat_flux=0. ! default value
207         call getin("tke_heat_flux",tke_heat_flux)
208         write(*,*) " tke_heat_flux = ",tke_heat_flux
209         write(*,*) " 0 means the usual schemes are computing"
210
211         write(*,*) "call radiative transfer ?"
212         callrad=.true. ! default value
213         call getin("callrad",callrad)
214         write(*,*) " callrad = ",callrad
215
216         write(*,*) "call slope insolation scheme ?",
217     &              "(matters only if callrad=T)"
218#ifdef MESOSCALE
219         callslope=.true. ! default value
220#else
221         callslope=.false. ! default value (not supported yet)
222#endif
223         call getin("callslope",callslope)
224         write(*,*) " callslope = ",callslope
225
226         write(*,*) "call NLTE radiative schemes ?",
227     &              "(matters only if callrad=T)"
228         callnlte=.false. ! default value
229         call getin("callnlte",callnlte)
230         write(*,*) " callnlte = ",callnlte
231         
232         nltemodel=0    !default value
233         write(*,*) "NLTE model?"
234         write(*,*) "0 -> old model, static O"
235         write(*,*) "1 -> old model, dynamic O"
236         write(*,*) "2 -> new model"
237         write(*,*) "(matters only if callnlte=T)"
238         call getin("nltemodel",nltemodel)
239         write(*,*) " nltemodel = ",nltemodel
240
241         write(*,*) "call CO2 NIR absorption ?",
242     &              "(matters only if callrad=T)"
243         callnirco2=.false. ! default value
244         call getin("callnirco2",callnirco2)
245         write(*,*) " callnirco2 = ",callnirco2
246
247         write(*,*) "New NIR NLTE correction ?",
248     $              "0-> old model (no correction)",
249     $              "1-> new correction",
250     $              "(matters only if callnirco2=T)"
251         nircorr=1      !default value
252         call getin("nircorr",nircorr)
253         write(*,*) " nircorr = ",nircorr
254
255         write(*,*) "call turbulent vertical diffusion ?"
256         calldifv=.true. ! default value
257         call getin("calldifv",calldifv)
258         write(*,*) " calldifv = ",calldifv
259
260         write(*,*) "call thermals ?"
261         calltherm=.false. ! default value
262         call getin("calltherm",calltherm)
263         write(*,*) " calltherm = ",calltherm
264
265         write(*,*) "output thermal diagnostics ?"
266         outptherm=.false. ! default value
267         call getin("outptherm",outptherm)
268         write(*,*) " outptherm = ",outptherm
269
270         write(*,*) "call convective adjustment ?"
271         calladj=.true. ! default value
272         call getin("calladj",calladj)
273         write(*,*) " calladj = ",calladj
274         
275         if (calltherm .and. (.not. calladj)) then
276          print*,'Convadj has to be activated when using thermals'
277          stop
278         endif
279
280         write(*,*) "call Richardson-based surface layer ?"
281         callrichsl=.false. ! default value
282         call getin("callrichsl",callrichsl)
283         write(*,*) " callrichsl = ",callrichsl
284
285         if (calltherm .and. .not.callrichsl) then
286          print*,'WARNING WARNING WARNING'
287          print*,'if calltherm=T we strongly advise that '
288          print*,'you use the new surface layer scheme '
289          print*,'by setting callrichsl=T '
290         endif
291
292         if (calladj .and. callrichsl .and. (.not. calltherm)) then
293          print*,'You should not be calling the convective adjustment
294     & scheme with the Richardson surface-layer and without the thermals
295     &. This approach is not
296     & physically consistent and can lead to unrealistic friction
297     & values.'
298          print*,'If you want to use the Ri. surface-layer, either
299     & activate thermals OR de-activate the convective adjustment.'
300          stop
301         endif
302
303         write(*,*) "call CO2 condensation ?"
304         callcond=.true. ! default value
305         call getin("callcond",callcond)
306         write(*,*) " callcond = ",callcond
307
308         write(*,*)"call thermal conduction in the soil ?"
309         callsoil=.true. ! default value
310         call getin("callsoil",callsoil)
311         write(*,*) " callsoil = ",callsoil
312         
313
314         write(*,*)"call Lott's gravity wave/subgrid topography ",
315     &             "scheme ?"
316         calllott=.true. ! default value
317         call getin("calllott",calllott)
318         write(*,*)" calllott = ",calllott
319
320
321         write(*,*)"rad.transfer is computed every iradia",
322     &             " physical timestep"
323         iradia=1 ! default value
324         call getin("iradia",iradia)
325         write(*,*)" iradia = ",iradia
326         
327
328         write(*,*)"Output of the exchange coefficient mattrix ?",
329     &             "(for diagnostics only)"
330         callg2d=.false. ! default value
331         call getin("callg2d",callg2d)
332         write(*,*)" callg2d = ",callg2d
333
334         write(*,*)"Rayleigh scattering : (should be .false. for now)"
335         rayleigh=.false.
336         call getin("rayleigh",rayleigh)
337         write(*,*)" rayleigh = ",rayleigh
338
339
340! TRACERS:
341
342! dustbin
343         write(*,*)"Transported dust ? (if >0, use 'dustbin' dust bins)"
344         dustbin=0 ! default value
345         call getin("dustbin",dustbin)
346         write(*,*)" dustbin = ",dustbin
347! active
348         write(*,*)"Radiatively active dust ? (matters if dustbin>0)"
349         active=.false. ! default value
350         call getin("active",active)
351         write(*,*)" active = ",active
352
353! Test of incompatibility:
354! if active is used, then dustbin should be > 0
355
356         if (active.and.(dustbin.lt.1)) then
357           print*,'if active is used, then dustbin should > 0'
358           stop
359         endif
360! doubleq
361         write(*,*)"use mass and number mixing ratios to predict",
362     &             " dust size ?"
363         doubleq=.false. ! default value
364         call getin("doubleq",doubleq)
365         write(*,*)" doubleq = ",doubleq
366! submicron
367         submicron=.false. ! default value
368         call getin("submicron",submicron)
369         write(*,*)" submicron = ",submicron
370
371! Test of incompatibility:
372! if doubleq is used, then dustbin should be 2
373
374         if (doubleq.and.(dustbin.ne.2)) then
375           print*,'if doubleq is used, then dustbin should be 2'
376           stop
377         endif
378         if (doubleq.and.submicron.and.(nqmx.LT.3)) then
379           print*,'If doubleq is used with a submicron tracer,'
380           print*,' then the number of tracers has to be'
381           print*,' larger than 3.'
382           stop
383         endif
384! lifting
385         write(*,*)"dust lifted by GCM surface winds ?"
386         lifting=.false. ! default value
387         call getin("lifting",lifting)
388         write(*,*)" lifting = ",lifting
389
390! Test of incompatibility:
391! if lifting is used, then dustbin should be > 0
392
393         if (lifting.and.(dustbin.lt.1)) then
394           print*,'if lifting is used, then dustbin should > 0'
395           stop
396         endif
397! callddevil
398         write(*,*)" dust lifted by dust devils ?"
399         callddevil=.false. !default value
400         call getin("callddevil",callddevil)
401         write(*,*)" callddevil = ",callddevil
402
403! Test of incompatibility:
404! if dustdevil is used, then dustbin should be > 0
405
406         if (callddevil.and.(dustbin.lt.1)) then
407           print*,'if dustdevil is used, then dustbin should > 0'
408           stop
409         endif
410! sedimentation
411         write(*,*) "Gravitationnal sedimentation ?"
412         sedimentation=.true. ! default value
413         call getin("sedimentation",sedimentation)
414         write(*,*) " sedimentation = ",sedimentation
415! activice
416         write(*,*) "Radiatively active transported atmospheric ",
417     &              "water ice ?"
418         activice=.false. ! default value
419         call getin("activice",activice)
420         write(*,*) " activice = ",activice
421! water
422         write(*,*) "Compute water cycle ?"
423         water=.false. ! default value
424         call getin("water",water)
425         write(*,*) " water = ",water
426
427! Test of incompatibility:
428
429         if (activice.and..not.water) then
430           print*,'if activice is used, water should be used too'
431           stop
432         endif
433
434         if (water.and..not.tracer) then
435           print*,'if water is used, tracer should be used too'
436           stop
437         endif
438         
439! water ice clouds effective variance distribution for sedimentaion       
440        write(*,*) "effective variance for water ice clouds ?"
441        nuice_sed=0.45
442        call getin("nuice_sed",nuice_sed)
443        write(*,*) "water_param nueff Sedimentation:", nuice_sed
444         
445! ccn factor if no scavenging         
446        write(*,*) "water param CCN reduc. factor ?", ccn_factor
447        ccn_factor = 4.5
448        call getin("ccn_factor",ccn_factor)
449        write(*,*)" ccn_factor = ",ccn_factor
450        write(*,*)"Careful: only used when microphys=F, otherwise"
451        write(*,*)"the contact parameter is used instead;"
452
453! microphys
454         write(*,*)"Microphysical scheme for water-ice clouds?"
455         microphys=.false. ! default value
456         call getin("microphys",microphys)
457         write(*,*)" microphys = ",microphys
458
459! microphysical parameter contact       
460         write(*,*) "water contact parameter ?"
461         mteta  = 0.95
462         call getin("mteta",mteta)
463         write(*,*) "mteta = ", mteta
464
465! scavenging
466         write(*,*)"Dust scavenging by H2O/CO2 snowfall ?"
467         scavenging=.false. ! default value
468         call getin("scavenging",scavenging)
469         write(*,*)" scavenging = ",scavenging
470         
471
472! Test of incompatibility:
473! if scavenging is used, then dustbin should be > 0
474
475         if (scavenging.and.(dustbin.lt.1)) then
476           print*,'if scavenging is used, then dustbin should > 0'
477           stop
478         endif
479         if ((microphys.and..not.doubleq).or.
480     &       (microphys.and..not.active).or.
481     &       (microphys.and..not.scavenging).or.
482     &       (microphys.and..not.water)) then
483           print*,'if microphys is used, then doubleq,'
484           print*,'active, water and scavenging must all be used!'
485           stop
486         endif
487         if ((scavenging.and..not.doubleq).or.
488     &       (scavenging.and..not.active).or.
489     &       (scavenging.and..not.water).or.
490     &       (scavenging.and..not.microphys)) then
491           print*,'if scavenging is used, then doubleq,'
492           print*,'active, water and microphys must all be used!'
493           stop
494         endif
495
496! Test of incompatibility:
497
498         write(*,*) "Permanent water caps at poles ?",
499     &               " .true. is RECOMMENDED"
500         write(*,*) "(with .true., North cap is a source of water ",
501     &   "and South pole is a cold trap)"
502         caps=.true. ! default value
503         call getin("caps",caps)
504         write(*,*) " caps = ",caps
505
506! albedo_h2o_ice
507         write(*,*) "water ice albedo ?"
508         albedo_h2o_ice=0.45
509         call getin("albedo_h2o_ice",albedo_h2o_ice)
510         write(*,*) " albedo_h2o_ice = ",albedo_h2o_ice
511! inert_h2o_ice
512         write(*,*) "water ice thermal inertia ?"
513         inert_h2o_ice=2400 ! (J.m^-2.K^-1.s^-1/2)
514         call getin("inert_h2o_ice",inert_h2o_ice)
515         write(*,*) " inert_h2o_ice = ",inert_h2o_ice
516! frost_albedo_threshold
517         write(*,*) "frost thickness threshold for albedo ?"
518         frost_albedo_threshold=0.005 ! 5.4 mic (i.e 0.005 kg.m-2)
519         call getin("frost_albedo_threshold",
520     &    frost_albedo_threshold)
521         write(*,*) " frost_albedo_threshold = ",
522     &            frost_albedo_threshold
523
524! call Titus crocus line -- DEFAULT IS NONE
525         write(*,*) "Titus crocus line ?"
526         tituscap=.false.  ! default value
527         call getin("tituscap",tituscap)
528         write(*,*) "tituscap",tituscap
529
530!!!!!!!!!!!!!!!! TEMP !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
531!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!     
532         write(*,*) "temp tag for water caps ?"
533         temptag=.false.
534         call getin("temptag",temptag)
535         write(*,*) " temptag = ",temptag
536!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
537!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!             
538
539         write(*,*) "photochemistry: include chemical species"
540         photochem=.false. ! default value
541         call getin("photochem",photochem)
542         write(*,*) " photochem = ",photochem
543
544
545! THERMOSPHERE
546
547         write(*,*) "call thermosphere ?"
548         callthermos=.false. ! default value
549         call getin("callthermos",callthermos)
550         write(*,*) " callthermos = ",callthermos
551         
552
553         write(*,*) " water included without cycle ",
554     &              "(only if water=.false.)"
555         thermoswater=.false. ! default value
556         call getin("thermoswater",thermoswater)
557         write(*,*) " thermoswater = ",thermoswater
558
559         write(*,*) "call thermal conduction ?",
560     &    " (only if callthermos=.true.)"
561         callconduct=.false. ! default value
562         call getin("callconduct",callconduct)
563         write(*,*) " callconduct = ",callconduct
564
565         write(*,*) "call EUV heating ?",
566     &   " (only if callthermos=.true.)"
567         calleuv=.false.  ! default value
568         call getin("calleuv",calleuv)
569         write(*,*) " calleuv = ",calleuv
570
571         write(*,*) "call molecular viscosity ?",
572     &   " (only if callthermos=.true.)"
573         callmolvis=.false. ! default value
574         call getin("callmolvis",callmolvis)
575         write(*,*) " callmolvis = ",callmolvis
576
577         write(*,*) "call molecular diffusion ?",
578     &   " (only if callthermos=.true.)"
579         callmoldiff=.false. ! default value
580         call getin("callmoldiff",callmoldiff)
581         write(*,*) " callmoldiff = ",callmoldiff
582         
583
584         write(*,*) "call thermospheric photochemistry ?",
585     &   " (only if callthermos=.true.)"
586         thermochem=.false. ! default value
587         call getin("thermochem",thermochem)
588         write(*,*) " thermochem = ",thermochem
589
590         write(*,*) "date for solar flux calculation:",
591     &   " (1985 < date < 2002)"
592         write(*,*) "(Solar min=1996.4 ave=1993.4 max=1990.6)"
593         solarcondate=1993.4 ! default value
594         call getin("solarcondate",solarcondate)
595         write(*,*) " solarcondate = ",solarcondate
596         
597         write(*,*) "UV heating efficiency:",
598     &   "measured values between 0.19 and 0.23 (Fox et al. 1996)",
599     &   "lower values may be used to compensate low 15 um cooling"
600         euveff=0.21 !default value
601         call getin("euveff",euveff)
602         write(*,*) " euveff = ", euveff
603
604         if (.not.callthermos) then
605           if (thermoswater) then
606             print*,'if thermoswater is set, callthermos must be true'
607             stop
608           endif         
609           if (callconduct) then
610             print*,'if callconduct is set, callthermos must be true'
611             stop
612           endif       
613           if (calleuv) then
614             print*,'if calleuv is set, callthermos must be true'
615             stop
616           endif         
617           if (callmolvis) then
618             print*,'if callmolvis is set, callthermos must be true'
619             stop
620           endif       
621           if (callmoldiff) then
622             print*,'if callmoldiff is set, callthermos must be true'
623             stop
624           endif         
625           if (thermochem) then
626             print*,'if thermochem is set, callthermos must be true'
627             stop
628           endif         
629        endif
630
631! Test of incompatibility:
632! if photochem is used, then water should be used too
633
634         if (photochem.and..not.water) then
635           print*,'if photochem is used, water should be used too'
636           stop
637         endif
638
639! if callthermos is used, then thermoswater should be used too
640! (if water not used already)
641
642         if (callthermos .and. .not.water) then
643           if (callthermos .and. .not.thermoswater) then
644             print*,'if callthermos is used, water or thermoswater
645     &               should be used too'
646             stop
647           endif
648         endif
649
650         PRINT*,'--------------------------------------------'
651         PRINT*
652         PRINT*
653      ELSE
654         write(*,*)
655         write(*,*) 'Cannot read file callphys.def. Is it here ?'
656         stop
657      ENDIF
658
6598000  FORMAT(t5,a12,l8)
6608001  FORMAT(t5,a12,i8)
661
662      PRINT*
663      PRINT*,'inifis: daysec',daysec
664      PRINT*
665      PRINT*,'inifis: The radiative transfer is computed:'
666      PRINT*,'           each ',iradia,' physical time-step'
667      PRINT*,'        or each ',iradia*dtphys,' seconds'
668      PRINT*
669! --------------------------------------------------------------
670!  Managing the Longwave radiative transfer
671! --------------------------------------------------------------
672
673!     In most cases, the run just use the following values :
674!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
675      callemis=.true.     
676!     ilwd=10*int(daysec/dtphys) ! bug before 22/10/01       
677      ilwd=1
678      ilwn=1 !2
679      ilwb=1 !2
680      linear=.true.       
681      ncouche=3
682      alphan=0.4
683      semi=0
684
685!     BUT people working hard on the LW may want to read them in 'radia.def'
686!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
687
688      OPEN(99,file='radia.def',status='old',form='formatted'
689     .     ,iostat=ierr)
690      IF(ierr.EQ.0) THEN
691         write(*,*) 'inifis: Reading radia.def !!!'
692         READ(99,fmt='(a)') ch1
693         READ(99,*) callemis
694         WRITE(*,8000) ch1,callemis
695
696         READ(99,fmt='(a)') ch1
697         READ(99,*) iradia
698         WRITE(*,8001) ch1,iradia
699
700         READ(99,fmt='(a)') ch1
701         READ(99,*) ilwd
702         WRITE(*,8001) ch1,ilwd
703
704         READ(99,fmt='(a)') ch1
705         READ(99,*) ilwn
706         WRITE(*,8001) ch1,ilwn
707
708         READ(99,fmt='(a)') ch1
709         READ(99,*) linear
710         WRITE(*,8000) ch1,linear
711
712         READ(99,fmt='(a)') ch1
713         READ(99,*) ncouche
714         WRITE(*,8001) ch1,ncouche
715
716         READ(99,fmt='(a)') ch1
717         READ(99,*) alphan
718         WRITE(*,*) ch1,alphan
719
720         READ(99,fmt='(a)') ch1
721         READ(99,*) ilwb
722         WRITE(*,8001) ch1,ilwb
723
724
725         READ(99,fmt='(a)') ch1
726         READ(99,'(l1)') callg2d
727         WRITE(*,8000) ch1,callg2d
728
729         READ(99,fmt='(a)') ch1
730         READ(99,*) semi
731         WRITE(*,*) ch1,semi
732      end if
733      CLOSE(99)
734
735!-----------------------------------------------------------------------
736!     Some more initialization:
737!     ------------------------
738
739      ! in 'comgeomfi.h'
740      CALL SCOPY(ngrid,plon,1,long,1)
741      CALL SCOPY(ngrid,plat,1,lati,1)
742      CALL SCOPY(ngrid,parea,1,area,1)
743      totarea=SSUM(ngridmx,area,1)
744
745      ! in 'comdiurn.h'
746      DO ig=1,ngrid
747         sinlat(ig)=sin(plat(ig))
748         coslat(ig)=cos(plat(ig))
749         sinlon(ig)=sin(plon(ig))
750         coslon(ig)=cos(plon(ig))
751      ENDDO
752
753      pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h
754
755!     managing the tracers, and tests:
756!     -------------------------------
757!     Ehouarn: removed; as these tests are now done in initracer.F
758!      if(tracer) then
759!
760!!          when photochem is used, nqchem_min is the rank
761!!          of the first chemical species
762!
763!! Ehouarn: nqchem_min is now meaningless and no longer used
764!!       nqchem_min = 1
765!       if (photochem .or. callthermos) then
766!         chem = .true.
767!       end if
768!
769!       if (water .or. thermoswater) h2o = .true.
770!
771!!          TESTS
772!
773!       print*,'inifis: TRACERS:'
774!       write(*,*) "    chem=",chem,"    h2o=",h2o
775!!       write(*,*) "   doubleq=",doubleq
776!!       write(*,*) "   dustbin=",dustbin
777!
778!       if ((doubleq).and.(h2o).and.
779!     $     (chem)) then
780!         print*,' 2 dust tracers (doubleq)'
781!         print*,' 1 water vapour tracer'
782!         print*,' 1 water ice tracer'
783!         print*,nqmx-4,' chemistry tracers'
784!       endif
785!
786!       if ((doubleq).and.(h2o).and.
787!     $     .not.(chem)) then
788!         print*,' 2 dust tracers (doubleq)'
789!         print*,' 1 water vapour tracer'
790!         print*,' 1 water ice tracer'
791!         if (nqmx.LT.4) then
792!           print*,'nqmx should be at least equal to'
793!           print*,'4 with these options.'
794!           stop
795!         endif
796!       endif
797!
798!       if (.not.(doubleq).and.(h2o).and.
799!     $     (chem)) then
800!         if (dustbin.gt.0) then
801!           print*,dustbin,' dust bins'
802!         endif
803!         print*,nqmx-2-dustbin,' chemistry tracers'
804!         print*,' 1 water vapour tracer'
805!         print*,' 1 water ice tracer'
806!       endif
807!
808!       if (.not.(doubleq).and.(h2o).and.
809!     $     .not.(chem)) then
810!         if (dustbin.gt.0) then
811!           print*,dustbin,' dust bins'
812!         endif
813!         print*,' 1 water vapour tracer'
814!         print*,' 1 water ice tracer'
815!         if (nqmx.gt.(dustbin+2)) then
816!           print*,'nqmx should be ',(dustbin+2),
817!     $            ' with these options...'
818!                  print*,'(or check callphys.def)'
819!         endif
820!         if (nqmx.lt.(dustbin+2)) then
821!           write(*,*) "inifis: nqmx.lt.(dustbin+2)"
822!           stop
823!         endif
824!       endif
825!
826!      endif ! of if (tracer)
827!
828!      RETURN
829      END
Note: See TracBrowser for help on using the repository browser.