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

Last change on this file since 308 was 288, checked in by aslmd, 13 years ago

LMDZ.MARS

09/09/11 == AS

--> Fixed a problem with .eq. used with booleans in physiq.F
[some good picky compilers complain about this]
--> Added a warning in inifis about using
callrichsl = .false.
calltherm = .true.
which is not recommended. The new surface layer has been built to go
with the new mixing layer scheme (thermals). And anyway it is a much
better scheme than the previous one.

File size: 23.1 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#ifdef MESOSCALE
64#include "comsoil.h"     !!MESOSCALE -- needed to fill volcapa
65#include "meso_inc/meso_inc_inifisvar.F"
66#endif
67      REAL prad,pg,pr,pcpp,pdaysec
68
69      REAL ptimestep
70      INTEGER day_ini
71
72      INTEGER ngrid,nlayer
73      REAL plat(ngrid),plon(ngrid),parea(ngridmx)
74      INTEGER ig,ierr
75 
76!      EXTERNAL iniorbit,orbite
77      EXTERNAL SSUM
78      REAL SSUM
79 
80      CHARACTER ch1*12
81      CHARACTER ch80*80
82
83!      logical chem, h2o
84
85!      chem = .false.
86!      h2o = .false.
87
88      rad=prad
89      cpp=pcpp
90      g=pg
91      r=pr
92      rcp=r/cpp
93      daysec=pdaysec
94      dtphys=ptimestep
95#ifdef MESOSCALE
96#include "meso_inc/meso_inc_inifisini.F"
97#endif
98
99! --------------------------------------------------------
100!     The usual Tests
101!     --------------
102      IF (nlayer.NE.nlayermx) THEN
103         PRINT*,'STOP in inifis'
104         PRINT*,'Probleme de dimensions :'
105         PRINT*,'nlayer     = ',nlayer
106         PRINT*,'nlayermx   = ',nlayermx
107         STOP
108      ENDIF
109
110      IF (ngrid.NE.ngridmx) THEN
111         PRINT*,'STOP in inifis'
112         PRINT*,'Probleme de dimensions :'
113         PRINT*,'ngrid     = ',ngrid
114         PRINT*,'ngridmx   = ',ngridmx
115         STOP
116      ENDIF
117
118! --------------------------------------------------------------
119!  Reading the "callphys.def" file controlling some key options
120! --------------------------------------------------------------
121     
122      ! check that 'callphys.def' file is around
123      OPEN(99,file='callphys.def',status='old',form='formatted'
124     &     ,iostat=ierr)
125      CLOSE(99)
126     
127      IF(ierr.EQ.0) THEN
128         PRINT*
129         PRINT*
130         PRINT*,'--------------------------------------------'
131         PRINT*,' inifis: Parameters for the physics (callphys.def)'
132         PRINT*,'--------------------------------------------'
133
134         write(*,*) "Directory where external input files are:"
135         datafile="/u/forget/WWW/datagcm/datafile"
136         call getin("datadir",datafile) ! default path
137         write(*,*) " datafile = ",trim(datafile)
138
139         write(*,*) "Run with or without tracer transport ?"
140         tracer=.false. ! default value
141         call getin("tracer",tracer)
142         write(*,*) " tracer = ",tracer
143
144         write(*,*) "Diurnal cycle ?"
145         write(*,*) "(if diurnal=False, diurnal averaged solar heating)"
146         diurnal=.true. ! default value
147         call getin("diurnal",diurnal)
148         write(*,*) " diurnal = ",diurnal
149
150         write(*,*) "Seasonal cycle ?"
151         write(*,*) "(if season=False, Ls stays constant, to value ",
152     &   "set in 'start'"
153         season=.true. ! default value
154         call getin("season",season)
155         write(*,*) " season = ",season
156
157         write(*,*) "Write some extra output to the screen ?"
158         lwrite=.false. ! default value
159         call getin("lwrite",lwrite)
160         write(*,*) " lwrite = ",lwrite
161
162         write(*,*) "Save statistics in file stats.nc ?"
163#ifdef MESOSCALE
164         callstats=.false. ! default value
165#else
166         callstats=.true. ! default value
167#endif
168         call getin("callstats",callstats)
169         write(*,*) " callstats = ",callstats
170
171         write(*,*) "Save EOF profiles in file 'profiles' for ",
172     &              "Climate Database?"
173         calleofdump=.false. ! default value
174         call getin("calleofdump",calleofdump)
175         write(*,*) " calleofdump = ",calleofdump
176
177         write(*,*) "Dust scenario: 1=constant dust (read from startfi",
178     &   " or set as tauvis); 2=Viking scenario; =3 MGS scenario,",
179     &   "=4 Mars Year 24 from TES assimilation, ",
180     &   "=24,25 or 26 :Mars Year 24,25 or 26 from TES assimilation"
181         iaervar=3 ! default value
182         call getin("iaervar",iaervar)
183         write(*,*) " iaervar = ",iaervar
184
185         write(*,*) "Reference (visible) dust opacity at 700 Pa ",
186     &   "(matters only if iaervar=1)"
187         ! NB: default value of tauvis is set/read in startfi.nc file
188         call getin("tauvis",tauvis)
189         write(*,*) " tauvis = ",tauvis
190
191         write(*,*) "Dust vertical distribution:"
192         write(*,*) "(=1 top set by topdustref parameter;",
193     & " =2 Viking scenario; =3 MGS scenario)"
194         iddist=3 ! default value
195         call getin("iddist",iddist)
196         write(*,*) " iddist = ",iddist
197
198         write(*,*) "Dust top altitude (km). (Matters only if iddist=1)"
199         topdustref= 90.0 ! default value
200         call getin("topdustref",topdustref)
201         write(*,*) " topdustref = ",topdustref
202
203         write(*,*) "call radiative transfer ?"
204         callrad=.true. ! default value
205         call getin("callrad",callrad)
206         write(*,*) " callrad = ",callrad
207
208         write(*,*) "call slope insolation scheme ?",
209     &              "(matters only if callrad=T)"
210#ifdef MESOSCALE
211         callslope=.true. ! default value
212#else
213         callslope=.false. ! default value (not supported yet)
214#endif
215         call getin("callslope",callslope)
216         write(*,*) " callslope = ",callslope
217
218         write(*,*) "call NLTE radiative schemes ?",
219     &              "(matters only if callrad=T)"
220         callnlte=.false. ! default value
221         call getin("callnlte",callnlte)
222         write(*,*) " callnlte = ",callnlte
223         
224         write(*,*) "call CO2 NIR absorption ?",
225     &              "(matters only if callrad=T)"
226         callnirco2=.false. ! default value
227         call getin("callnirco2",callnirco2)
228         write(*,*) " callnirco2 = ",callnirco2
229         
230         write(*,*) "call turbulent vertical diffusion ?"
231         calldifv=.true. ! default value
232         call getin("calldifv",calldifv)
233         write(*,*) " calldifv = ",calldifv
234
235         write(*,*) "call thermals ?"
236         calltherm=.false. ! default value
237         call getin("calltherm",calltherm)
238         write(*,*) " calltherm = ",calltherm
239
240         write(*,*) "output thermal diagnostics ?"
241         outptherm=.false. ! default value
242         call getin("outptherm",outptherm)
243         write(*,*) " outptherm = ",outptherm
244
245         write(*,*) "call convective adjustment ?"
246         calladj=.true. ! default value
247         call getin("calladj",calladj)
248         write(*,*) " calladj = ",calladj
249         
250         if (calltherm .and. (.not. calladj)) then
251          print*,'Convadj has to be activated when using thermals'
252          stop
253         endif
254
255         write(*,*) "call Richardson-based surface layer ?"
256         callrichsl=.false. ! default value
257         call getin("callrichsl",callrichsl)
258         write(*,*) " callrichsl = ",callrichsl
259
260         if (calltherm .and. .not.callrichsl) then
261          print*,'WARNING WARNING WARNING'
262          print*,'if calltherm=T we strongly advise that '
263          print*,'you use the new surface layer scheme '
264          print*,'by setting callrichsl=T '
265         endif
266
267         write(*,*) "call CO2 condensation ?"
268         callcond=.true. ! default value
269         call getin("callcond",callcond)
270         write(*,*) " callcond = ",callcond
271
272         write(*,*)"call thermal conduction in the soil ?"
273         callsoil=.true. ! default value
274         call getin("callsoil",callsoil)
275         write(*,*) " callsoil = ",callsoil
276         
277
278         write(*,*)"call Lott's gravity wave/subgrid topography ",
279     &             "scheme ?"
280         calllott=.true. ! default value
281         call getin("calllott",calllott)
282         write(*,*)" calllott = ",calllott
283
284
285         write(*,*)"rad.transfer is computed every iradia",
286     &             " physical timestep"
287         iradia=1 ! default value
288         call getin("iradia",iradia)
289         write(*,*)" iradia = ",iradia
290         
291
292         write(*,*)"Output of the exchange coefficient mattrix ?",
293     &             "(for diagnostics only)"
294         callg2d=.false. ! default value
295         call getin("callg2d",callg2d)
296         write(*,*)" callg2d = ",callg2d
297
298         write(*,*)"Rayleigh scattering : (should be .false. for now)"
299         rayleigh=.false.
300         call getin("rayleigh",rayleigh)
301         write(*,*)" rayleigh = ",rayleigh
302
303
304! TRACERS:
305
306         write(*,*)"Transported dust ? (if >0, use 'dustbin' dust bins)"
307         dustbin=0 ! default value
308         call getin("dustbin",dustbin)
309         write(*,*)" dustbin = ",dustbin
310
311         write(*,*)"Radiatively active dust ? (matters if dustbin>0)"
312         active=.false. ! default value
313         call getin("active",active)
314         write(*,*)" active = ",active
315
316! Test of incompatibility:
317! if active is used, then dustbin should be > 0
318
319         if (active.and.(dustbin.lt.1)) then
320           print*,'if active is used, then dustbin should > 0'
321           stop
322         endif
323
324         write(*,*)"use mass and number mixing ratios to predict",
325     &             " dust size ?"
326         doubleq=.false. ! default value
327         call getin("doubleq",doubleq)
328         write(*,*)" doubleq = ",doubleq
329
330         submicron=.false. ! default value
331         call getin("submicron",submicron)
332         write(*,*)" submicron = ",submicron
333
334! Test of incompatibility:
335! if doubleq is used, then dustbin should be 2
336
337         if (doubleq.and.(dustbin.ne.2)) then
338           print*,'if doubleq is used, then dustbin should be 2'
339           stop
340         endif
341         if (doubleq.and.submicron.and.(nqmx.LT.3)) then
342           print*,'If doubleq is used with a submicron tracer,'
343           print*,' then the number of tracers has to be'
344           print*,' larger than 3.'
345           stop
346         endif
347
348         write(*,*)"dust lifted by GCM surface winds ?"
349         lifting=.false. ! default value
350         call getin("lifting",lifting)
351         write(*,*)" lifting = ",lifting
352
353! Test of incompatibility:
354! if lifting is used, then dustbin should be > 0
355
356         if (lifting.and.(dustbin.lt.1)) then
357           print*,'if lifting is used, then dustbin should > 0'
358           stop
359         endif
360
361         write(*,*)" dust lifted by dust devils ?"
362         callddevil=.false. !default value
363         call getin("callddevil",callddevil)
364         write(*,*)" callddevil = ",callddevil
365         
366
367! Test of incompatibility:
368! if dustdevil is used, then dustbin should be > 0
369
370         if (callddevil.and.(dustbin.lt.1)) then
371           print*,'if dustdevil is used, then dustbin should > 0'
372           stop
373         endif
374
375         write(*,*)"Dust scavenging by CO2 snowfall ?"
376         scavenging=.false. ! default value
377         call getin("scavenging",scavenging)
378         write(*,*)" scavenging = ",scavenging
379         
380
381! Test of incompatibility:
382! if scavenging is used, then dustbin should be > 0
383
384         if (scavenging.and.(dustbin.lt.1)) then
385           print*,'if scavenging is used, then dustbin should > 0'
386           stop
387         endif
388
389         write(*,*) "Gravitationnal sedimentation ?"
390         sedimentation=.true. ! default value
391         call getin("sedimentation",sedimentation)
392         write(*,*) " sedimentation = ",sedimentation
393
394         write(*,*) "Radiatively active transported atmospheric ",
395     &              "water ice ?"
396         activice=.false. ! default value
397         call getin("activice",activice)
398         write(*,*) " activice = ",activice
399
400         write(*,*) "Compute water cycle ?"
401         water=.false. ! default value
402         call getin("water",water)
403         write(*,*) " water = ",water
404
405! Test of incompatibility:
406
407         if (activice.and..not.water) then
408           print*,'if activice is used, water should be used too'
409           stop
410         endif
411
412         if (water.and..not.tracer) then
413           print*,'if water is used, tracer should be used too'
414           stop
415         endif
416
417! Test of incompatibility:
418
419         write(*,*) "Permanent water caps at poles ?",
420     &               " .true. is RECOMMENDED"
421         write(*,*) "(with .true., North cap is a source of water ",
422     &   "and South pole is a cold trap)"
423         caps=.true. ! default value
424         call getin("caps",caps)
425         write(*,*) " caps = ",caps
426
427
428         write(*,*) "water ice albedo ?"
429         albedo_h2o_ice=0.45
430         call getin("albedo_h2o_ice",albedo_h2o_ice)
431         write(*,*) " albedo_h2o_ice = ",albedo_h2o_ice
432         
433         write(*,*) "water ice thermal inertia ?"
434         inert_h2o_ice=2400 ! (J.m^-2.K^-1.s^-1/2)
435         call getin("inert_h2o_ice",inert_h2o_ice)
436         write(*,*) " inert_h2o_ice = ",inert_h2o_ice
437         
438         write(*,*) "frost thickness threshold for albedo ?"
439         frost_albedo_threshold=0.005 ! 5.4 mic (i.e 0.005 kg.m-2)
440         call getin("frost_albedo_threshold",
441     &    frost_albedo_threshold)
442         write(*,*) " frost_albedo_threshold = ",
443     &            frost_albedo_threshold
444
445!!!!!!!!!!!!!!!! TEMP !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
446!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!     
447         write(*,*) "temp tag for water caps ?"
448         temptag=.false.
449         call getin("temptag",temptag)
450         write(*,*) " temptag = ",temptag
451!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
452!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!             
453
454         write(*,*) "photochemistry: include chemical species"
455         photochem=.false. ! default value
456         call getin("photochem",photochem)
457         write(*,*) " photochem = ",photochem
458
459
460! THERMOSPHERE
461
462         write(*,*) "call thermosphere ?"
463         callthermos=.false. ! default value
464         call getin("callthermos",callthermos)
465         write(*,*) " callthermos = ",callthermos
466         
467
468         write(*,*) " water included without cycle ",
469     &              "(only if water=.false.)"
470         thermoswater=.false. ! default value
471         call getin("thermoswater",thermoswater)
472         write(*,*) " thermoswater = ",thermoswater
473
474         write(*,*) "call thermal conduction ?",
475     &    " (only if callthermos=.true.)"
476         callconduct=.false. ! default value
477         call getin("callconduct",callconduct)
478         write(*,*) " callconduct = ",callconduct
479
480         write(*,*) "call EUV heating ?",
481     &   " (only if callthermos=.true.)"
482         calleuv=.false.  ! default value
483         call getin("calleuv",calleuv)
484         write(*,*) " calleuv = ",calleuv
485
486         write(*,*) "call molecular viscosity ?",
487     &   " (only if callthermos=.true.)"
488         callmolvis=.false. ! default value
489         call getin("callmolvis",callmolvis)
490         write(*,*) " callmolvis = ",callmolvis
491
492         write(*,*) "call molecular diffusion ?",
493     &   " (only if callthermos=.true.)"
494         callmoldiff=.false. ! default value
495         call getin("callmoldiff",callmoldiff)
496         write(*,*) " callmoldiff = ",callmoldiff
497         
498
499         write(*,*) "call thermospheric photochemistry ?",
500     &   " (only if callthermos=.true.)"
501         thermochem=.false. ! default value
502         call getin("thermochem",thermochem)
503         write(*,*) " thermochem = ",thermochem
504
505         write(*,*) "date for solar flux calculation:",
506     &   " (1985 < date < 2002)"
507         write(*,*) "(Solar min=1996.4 ave=1993.4 max=1990.6)"
508         solarcondate=1993.4 ! default value
509         call getin("solarcondate",solarcondate)
510         write(*,*) " solarcondate = ",solarcondate
511         
512
513         if (.not.callthermos) then
514           if (thermoswater) then
515             print*,'if thermoswater is set, callthermos must be true'
516             stop
517           endif         
518           if (callconduct) then
519             print*,'if callconduct is set, callthermos must be true'
520             stop
521           endif       
522           if (calleuv) then
523             print*,'if calleuv is set, callthermos must be true'
524             stop
525           endif         
526           if (callmolvis) then
527             print*,'if callmolvis is set, callthermos must be true'
528             stop
529           endif       
530           if (callmoldiff) then
531             print*,'if callmoldiff is set, callthermos must be true'
532             stop
533           endif         
534           if (thermochem) then
535             print*,'if thermochem is set, callthermos must be true'
536             stop
537           endif         
538        endif
539
540! Test of incompatibility:
541! if photochem is used, then water should be used too
542
543         if (photochem.and..not.water) then
544           print*,'if photochem is used, water should be used too'
545           stop
546         endif
547
548! if callthermos is used, then thermoswater should be used too
549! (if water not used already)
550
551         if (callthermos .and. .not.water) then
552           if (callthermos .and. .not.thermoswater) then
553             print*,'if callthermos is used, water or thermoswater
554     &               should be used too'
555             stop
556           endif
557         endif
558
559         PRINT*,'--------------------------------------------'
560         PRINT*
561         PRINT*
562      ELSE
563         write(*,*)
564         write(*,*) 'Cannot read file callphys.def. Is it here ?'
565         stop
566      ENDIF
567
5688000  FORMAT(t5,a12,l8)
5698001  FORMAT(t5,a12,i8)
570
571      PRINT*
572      PRINT*,'inifis: daysec',daysec
573      PRINT*
574      PRINT*,'inifis: The radiative transfer is computed:'
575      PRINT*,'           each ',iradia,' physical time-step'
576      PRINT*,'        or each ',iradia*dtphys,' seconds'
577      PRINT*
578! --------------------------------------------------------------
579!  Managing the Longwave radiative transfer
580! --------------------------------------------------------------
581
582!     In most cases, the run just use the following values :
583!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
584      callemis=.true.     
585!     ilwd=10*int(daysec/dtphys) ! bug before 22/10/01       
586      ilwd=1
587      ilwn=1 !2
588      ilwb=1 !2
589      linear=.true.       
590      ncouche=3
591      alphan=0.4
592      semi=0
593
594!     BUT people working hard on the LW may want to read them in 'radia.def'
595!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
596
597      OPEN(99,file='radia.def',status='old',form='formatted'
598     .     ,iostat=ierr)
599      IF(ierr.EQ.0) THEN
600         write(*,*) 'inifis: Reading radia.def !!!'
601         READ(99,fmt='(a)') ch1
602         READ(99,*) callemis
603         WRITE(*,8000) ch1,callemis
604
605         READ(99,fmt='(a)') ch1
606         READ(99,*) iradia
607         WRITE(*,8001) ch1,iradia
608
609         READ(99,fmt='(a)') ch1
610         READ(99,*) ilwd
611         WRITE(*,8001) ch1,ilwd
612
613         READ(99,fmt='(a)') ch1
614         READ(99,*) ilwn
615         WRITE(*,8001) ch1,ilwn
616
617         READ(99,fmt='(a)') ch1
618         READ(99,*) linear
619         WRITE(*,8000) ch1,linear
620
621         READ(99,fmt='(a)') ch1
622         READ(99,*) ncouche
623         WRITE(*,8001) ch1,ncouche
624
625         READ(99,fmt='(a)') ch1
626         READ(99,*) alphan
627         WRITE(*,*) ch1,alphan
628
629         READ(99,fmt='(a)') ch1
630         READ(99,*) ilwb
631         WRITE(*,8001) ch1,ilwb
632
633
634         READ(99,fmt='(a)') ch1
635         READ(99,'(l1)') callg2d
636         WRITE(*,8000) ch1,callg2d
637
638         READ(99,fmt='(a)') ch1
639         READ(99,*) semi
640         WRITE(*,*) ch1,semi
641      end if
642      CLOSE(99)
643
644!-----------------------------------------------------------------------
645!     Some more initialization:
646!     ------------------------
647
648      ! in 'comgeomfi.h'
649      CALL SCOPY(ngrid,plon,1,long,1)
650      CALL SCOPY(ngrid,plat,1,lati,1)
651      CALL SCOPY(ngrid,parea,1,area,1)
652      totarea=SSUM(ngridmx,area,1)
653
654      ! in 'comdiurn.h'
655      DO ig=1,ngrid
656         sinlat(ig)=sin(plat(ig))
657         coslat(ig)=cos(plat(ig))
658         sinlon(ig)=sin(plon(ig))
659         coslon(ig)=cos(plon(ig))
660      ENDDO
661
662      pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h
663
664!     managing the tracers, and tests:
665!     -------------------------------
666!     Ehouarn: removed; as these tests are now done in initracer.F
667!      if(tracer) then
668!
669!!          when photochem is used, nqchem_min is the rank
670!!          of the first chemical species
671!
672!! Ehouarn: nqchem_min is now meaningless and no longer used
673!!       nqchem_min = 1
674!       if (photochem .or. callthermos) then
675!         chem = .true.
676!       end if
677!
678!       if (water .or. thermoswater) h2o = .true.
679!
680!!          TESTS
681!
682!       print*,'inifis: TRACERS:'
683!       write(*,*) "    chem=",chem,"    h2o=",h2o
684!!       write(*,*) "   doubleq=",doubleq
685!!       write(*,*) "   dustbin=",dustbin
686!
687!       if ((doubleq).and.(h2o).and.
688!     $     (chem)) then
689!         print*,' 2 dust tracers (doubleq)'
690!         print*,' 1 water vapour tracer'
691!         print*,' 1 water ice tracer'
692!         print*,nqmx-4,' chemistry tracers'
693!       endif
694!
695!       if ((doubleq).and.(h2o).and.
696!     $     .not.(chem)) then
697!         print*,' 2 dust tracers (doubleq)'
698!         print*,' 1 water vapour tracer'
699!         print*,' 1 water ice tracer'
700!         if (nqmx.LT.4) then
701!           print*,'nqmx should be at least equal to'
702!           print*,'4 with these options.'
703!           stop
704!         endif
705!       endif
706!
707!       if (.not.(doubleq).and.(h2o).and.
708!     $     (chem)) then
709!         if (dustbin.gt.0) then
710!           print*,dustbin,' dust bins'
711!         endif
712!         print*,nqmx-2-dustbin,' chemistry tracers'
713!         print*,' 1 water vapour tracer'
714!         print*,' 1 water ice tracer'
715!       endif
716!
717!       if (.not.(doubleq).and.(h2o).and.
718!     $     .not.(chem)) then
719!         if (dustbin.gt.0) then
720!           print*,dustbin,' dust bins'
721!         endif
722!         print*,' 1 water vapour tracer'
723!         print*,' 1 water ice tracer'
724!         if (nqmx.gt.(dustbin+2)) then
725!           print*,'nqmx should be ',(dustbin+2),
726!     $            ' with these options...'
727!                  print*,'(or check callphys.def)'
728!         endif
729!         if (nqmx.lt.(dustbin+2)) then
730!           write(*,*) "inifis: nqmx.lt.(dustbin+2)"
731!           stop
732!         endif
733!       endif
734!
735!      endif ! of if (tracer)
736!
737!      RETURN
738      END
Note: See TracBrowser for help on using the repository browser.