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

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

LMDZ.MARS:

AS: J'ai fait ce commit et fait a la main un merge avec les modifications entre temps

24/08/11 == TN

Attempts to tune the water cycle by adding outliers

+ A few structural changes !!

  • watercap.h is now obsolete and removed -- all is in surfdat.h
  • watercaptag initialized in surfini.F (up to 7 areas defined) instead of initracer.F
    • settings proposed by AS commented
    • experiments by TN decommented. use with caution.
  • water ice albedo and thermal inertia in callphys.def and inifis.F
  • water ice albedo in surfini.F
  • water ice albedo computation in albedocaps.F90
  • alb_surfice is now obsolete in physiq.F, albedo_h2o_ice is used instead
  • frost_albedo_threshold defined in surfdat.h
  • water ice thermal inertia in soil.F

TODO: * calibrate thermal inertia and ice albedo

  • have a look at subgrid-scale ice with dryness ?
File size: 22.7 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 CO2 condensation ?"
256         callcond=.true. ! default value
257         call getin("callcond",callcond)
258         write(*,*) " callcond = ",callcond
259
260         write(*,*)"call thermal conduction in the soil ?"
261         callsoil=.true. ! default value
262         call getin("callsoil",callsoil)
263         write(*,*) " callsoil = ",callsoil
264         
265
266         write(*,*)"call Lott's gravity wave/subgrid topography ",
267     &             "scheme ?"
268         calllott=.true. ! default value
269         call getin("calllott",calllott)
270         write(*,*)" calllott = ",calllott
271
272
273         write(*,*)"rad.transfer is computed every iradia",
274     &             " physical timestep"
275         iradia=1 ! default value
276         call getin("iradia",iradia)
277         write(*,*)" iradia = ",iradia
278         
279
280         write(*,*)"Output of the exchange coefficient mattrix ?",
281     &             "(for diagnostics only)"
282         callg2d=.false. ! default value
283         call getin("callg2d",callg2d)
284         write(*,*)" callg2d = ",callg2d
285
286         write(*,*)"Rayleigh scattering : (should be .false. for now)"
287         rayleigh=.false.
288         call getin("rayleigh",rayleigh)
289         write(*,*)" rayleigh = ",rayleigh
290
291
292! TRACERS:
293
294         write(*,*)"Transported dust ? (if >0, use 'dustbin' dust bins)"
295         dustbin=0 ! default value
296         call getin("dustbin",dustbin)
297         write(*,*)" dustbin = ",dustbin
298
299         write(*,*)"Radiatively active dust ? (matters if dustbin>0)"
300         active=.false. ! default value
301         call getin("active",active)
302         write(*,*)" active = ",active
303
304! Test of incompatibility:
305! if active is used, then dustbin should be > 0
306
307         if (active.and.(dustbin.lt.1)) then
308           print*,'if active is used, then dustbin should > 0'
309           stop
310         endif
311
312         write(*,*)"use mass and number mixing ratios to predict",
313     &             " dust size ?"
314         doubleq=.false. ! default value
315         call getin("doubleq",doubleq)
316         write(*,*)" doubleq = ",doubleq
317
318         submicron=.false. ! default value
319         call getin("submicron",submicron)
320         write(*,*)" submicron = ",submicron
321
322! Test of incompatibility:
323! if doubleq is used, then dustbin should be 2
324
325         if (doubleq.and.(dustbin.ne.2)) then
326           print*,'if doubleq is used, then dustbin should be 2'
327           stop
328         endif
329         if (doubleq.and.submicron.and.(nqmx.LT.3)) then
330           print*,'If doubleq is used with a submicron tracer,'
331           print*,' then the number of tracers has to be'
332           print*,' larger than 3.'
333           stop
334         endif
335
336         write(*,*)"dust lifted by GCM surface winds ?"
337         lifting=.false. ! default value
338         call getin("lifting",lifting)
339         write(*,*)" lifting = ",lifting
340
341! Test of incompatibility:
342! if lifting is used, then dustbin should be > 0
343
344         if (lifting.and.(dustbin.lt.1)) then
345           print*,'if lifting is used, then dustbin should > 0'
346           stop
347         endif
348
349         write(*,*)" dust lifted by dust devils ?"
350         callddevil=.false. !default value
351         call getin("callddevil",callddevil)
352         write(*,*)" callddevil = ",callddevil
353         
354
355! Test of incompatibility:
356! if dustdevil is used, then dustbin should be > 0
357
358         if (callddevil.and.(dustbin.lt.1)) then
359           print*,'if dustdevil is used, then dustbin should > 0'
360           stop
361         endif
362
363         write(*,*)"Dust scavenging by CO2 snowfall ?"
364         scavenging=.false. ! default value
365         call getin("scavenging",scavenging)
366         write(*,*)" scavenging = ",scavenging
367         
368
369! Test of incompatibility:
370! if scavenging is used, then dustbin should be > 0
371
372         if (scavenging.and.(dustbin.lt.1)) then
373           print*,'if scavenging is used, then dustbin should > 0'
374           stop
375         endif
376
377         write(*,*) "Gravitationnal sedimentation ?"
378         sedimentation=.true. ! default value
379         call getin("sedimentation",sedimentation)
380         write(*,*) " sedimentation = ",sedimentation
381
382         write(*,*) "Radiatively active transported atmospheric ",
383     &              "water ice ?"
384         activice=.false. ! default value
385         call getin("activice",activice)
386         write(*,*) " activice = ",activice
387
388         write(*,*) "Compute water cycle ?"
389         water=.false. ! default value
390         call getin("water",water)
391         write(*,*) " water = ",water
392
393! Test of incompatibility:
394
395         if (activice.and..not.water) then
396           print*,'if activice is used, water should be used too'
397           stop
398         endif
399
400         if (water.and..not.tracer) then
401           print*,'if water is used, tracer should be used too'
402           stop
403         endif
404
405! Test of incompatibility:
406
407         write(*,*) "Permanent water caps at poles ?",
408     &               " .true. is RECOMMENDED"
409         write(*,*) "(with .true., North cap is a source of water ",
410     &   "and South pole is a cold trap)"
411         caps=.true. ! default value
412         call getin("caps",caps)
413         write(*,*) " caps = ",caps
414
415
416         write(*,*) "water ice albedo ?"
417         albedo_h2o_ice=0.45
418         call getin("albedo_h2o_ice",albedo_h2o_ice)
419         write(*,*) " albedo_h2o_ice = ",albedo_h2o_ice
420         
421         write(*,*) "water ice thermal inertia ?"
422         inert_h2o_ice=2400 ! (J.m^-2.K^-1.s^-1/2)
423         call getin("inert_h2o_ice",inert_h2o_ice)
424         write(*,*) " inert_h2o_ice = ",inert_h2o_ice
425         
426         write(*,*) "frost thickness threshold for albedo ?"
427         frost_albedo_threshold=0.005 ! 5.4 µm (i.e 0.005 kg/m²)
428         call getin("frost_albedo_threshold",
429     &    frost_albedo_threshold)
430         write(*,*) " frost_albedo_threshold = ",
431     &            frost_albedo_threshold
432
433!!!!!!!!!!!!!!!! TEMP !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
434!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!     
435         write(*,*) "temp tag for water caps ?"
436         temptag=.false.
437         call getin("temptag",temptag)
438         write(*,*) " temptag = ",temptag
439!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
440!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!             
441
442         write(*,*) "photochemistry: include chemical species"
443         photochem=.false. ! default value
444         call getin("photochem",photochem)
445         write(*,*) " photochem = ",photochem
446
447
448! THERMOSPHERE
449
450         write(*,*) "call thermosphere ?"
451         callthermos=.false. ! default value
452         call getin("callthermos",callthermos)
453         write(*,*) " callthermos = ",callthermos
454         
455
456         write(*,*) " water included without cycle ",
457     &              "(only if water=.false.)"
458         thermoswater=.false. ! default value
459         call getin("thermoswater",thermoswater)
460         write(*,*) " thermoswater = ",thermoswater
461
462         write(*,*) "call thermal conduction ?",
463     &    " (only if callthermos=.true.)"
464         callconduct=.false. ! default value
465         call getin("callconduct",callconduct)
466         write(*,*) " callconduct = ",callconduct
467
468         write(*,*) "call EUV heating ?",
469     &   " (only if callthermos=.true.)"
470         calleuv=.false.  ! default value
471         call getin("calleuv",calleuv)
472         write(*,*) " calleuv = ",calleuv
473
474         write(*,*) "call molecular viscosity ?",
475     &   " (only if callthermos=.true.)"
476         callmolvis=.false. ! default value
477         call getin("callmolvis",callmolvis)
478         write(*,*) " callmolvis = ",callmolvis
479
480         write(*,*) "call molecular diffusion ?",
481     &   " (only if callthermos=.true.)"
482         callmoldiff=.false. ! default value
483         call getin("callmoldiff",callmoldiff)
484         write(*,*) " callmoldiff = ",callmoldiff
485         
486
487         write(*,*) "call thermospheric photochemistry ?",
488     &   " (only if callthermos=.true.)"
489         thermochem=.false. ! default value
490         call getin("thermochem",thermochem)
491         write(*,*) " thermochem = ",thermochem
492
493         write(*,*) "date for solar flux calculation:",
494     &   " (1985 < date < 2002)"
495         write(*,*) "(Solar min=1996.4 ave=1993.4 max=1990.6)"
496         solarcondate=1993.4 ! default value
497         call getin("solarcondate",solarcondate)
498         write(*,*) " solarcondate = ",solarcondate
499         
500
501         if (.not.callthermos) then
502           if (thermoswater) then
503             print*,'if thermoswater is set, callthermos must be true'
504             stop
505           endif         
506           if (callconduct) then
507             print*,'if callconduct is set, callthermos must be true'
508             stop
509           endif       
510           if (calleuv) then
511             print*,'if calleuv is set, callthermos must be true'
512             stop
513           endif         
514           if (callmolvis) then
515             print*,'if callmolvis is set, callthermos must be true'
516             stop
517           endif       
518           if (callmoldiff) then
519             print*,'if callmoldiff is set, callthermos must be true'
520             stop
521           endif         
522           if (thermochem) then
523             print*,'if thermochem is set, callthermos must be true'
524             stop
525           endif         
526        endif
527
528! Test of incompatibility:
529! if photochem is used, then water should be used too
530
531         if (photochem.and..not.water) then
532           print*,'if photochem is used, water should be used too'
533           stop
534         endif
535
536! if callthermos is used, then thermoswater should be used too
537! (if water not used already)
538
539         if (callthermos .and. .not.water) then
540           if (callthermos .and. .not.thermoswater) then
541             print*,'if callthermos is used, water or thermoswater
542     &               should be used too'
543             stop
544           endif
545         endif
546
547         PRINT*,'--------------------------------------------'
548         PRINT*
549         PRINT*
550      ELSE
551         write(*,*)
552         write(*,*) 'Cannot read file callphys.def. Is it here ?'
553         stop
554      ENDIF
555
5568000  FORMAT(t5,a12,l8)
5578001  FORMAT(t5,a12,i8)
558
559      PRINT*
560      PRINT*,'inifis: daysec',daysec
561      PRINT*
562      PRINT*,'inifis: The radiative transfer is computed:'
563      PRINT*,'           each ',iradia,' physical time-step'
564      PRINT*,'        or each ',iradia*dtphys,' seconds'
565      PRINT*
566! --------------------------------------------------------------
567!  Managing the Longwave radiative transfer
568! --------------------------------------------------------------
569
570!     In most cases, the run just use the following values :
571!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
572      callemis=.true.     
573!     ilwd=10*int(daysec/dtphys) ! bug before 22/10/01       
574      ilwd=1
575      ilwn=1 !2
576      ilwb=1 !2
577      linear=.true.       
578      ncouche=3
579      alphan=0.4
580      semi=0
581
582!     BUT people working hard on the LW may want to read them in 'radia.def'
583!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
584
585      OPEN(99,file='radia.def',status='old',form='formatted'
586     .     ,iostat=ierr)
587      IF(ierr.EQ.0) THEN
588         write(*,*) 'inifis: Reading radia.def !!!'
589         READ(99,fmt='(a)') ch1
590         READ(99,*) callemis
591         WRITE(*,8000) ch1,callemis
592
593         READ(99,fmt='(a)') ch1
594         READ(99,*) iradia
595         WRITE(*,8001) ch1,iradia
596
597         READ(99,fmt='(a)') ch1
598         READ(99,*) ilwd
599         WRITE(*,8001) ch1,ilwd
600
601         READ(99,fmt='(a)') ch1
602         READ(99,*) ilwn
603         WRITE(*,8001) ch1,ilwn
604
605         READ(99,fmt='(a)') ch1
606         READ(99,*) linear
607         WRITE(*,8000) ch1,linear
608
609         READ(99,fmt='(a)') ch1
610         READ(99,*) ncouche
611         WRITE(*,8001) ch1,ncouche
612
613         READ(99,fmt='(a)') ch1
614         READ(99,*) alphan
615         WRITE(*,*) ch1,alphan
616
617         READ(99,fmt='(a)') ch1
618         READ(99,*) ilwb
619         WRITE(*,8001) ch1,ilwb
620
621
622         READ(99,fmt='(a)') ch1
623         READ(99,'(l1)') callg2d
624         WRITE(*,8000) ch1,callg2d
625
626         READ(99,fmt='(a)') ch1
627         READ(99,*) semi
628         WRITE(*,*) ch1,semi
629      end if
630      CLOSE(99)
631
632!-----------------------------------------------------------------------
633!     Some more initialization:
634!     ------------------------
635
636      ! in 'comgeomfi.h'
637      CALL SCOPY(ngrid,plon,1,long,1)
638      CALL SCOPY(ngrid,plat,1,lati,1)
639      CALL SCOPY(ngrid,parea,1,area,1)
640      totarea=SSUM(ngridmx,area,1)
641
642      ! in 'comdiurn.h'
643      DO ig=1,ngrid
644         sinlat(ig)=sin(plat(ig))
645         coslat(ig)=cos(plat(ig))
646         sinlon(ig)=sin(plon(ig))
647         coslon(ig)=cos(plon(ig))
648      ENDDO
649
650      pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h
651
652!     managing the tracers, and tests:
653!     -------------------------------
654!     Ehouarn: removed; as these tests are now done in initracer.F
655!      if(tracer) then
656!
657!!          when photochem is used, nqchem_min is the rank
658!!          of the first chemical species
659!
660!! Ehouarn: nqchem_min is now meaningless and no longer used
661!!       nqchem_min = 1
662!       if (photochem .or. callthermos) then
663!         chem = .true.
664!       end if
665!
666!       if (water .or. thermoswater) h2o = .true.
667!
668!!          TESTS
669!
670!       print*,'inifis: TRACERS:'
671!       write(*,*) "    chem=",chem,"    h2o=",h2o
672!!       write(*,*) "   doubleq=",doubleq
673!!       write(*,*) "   dustbin=",dustbin
674!
675!       if ((doubleq).and.(h2o).and.
676!     $     (chem)) then
677!         print*,' 2 dust tracers (doubleq)'
678!         print*,' 1 water vapour tracer'
679!         print*,' 1 water ice tracer'
680!         print*,nqmx-4,' chemistry tracers'
681!       endif
682!
683!       if ((doubleq).and.(h2o).and.
684!     $     .not.(chem)) then
685!         print*,' 2 dust tracers (doubleq)'
686!         print*,' 1 water vapour tracer'
687!         print*,' 1 water ice tracer'
688!         if (nqmx.LT.4) then
689!           print*,'nqmx should be at least equal to'
690!           print*,'4 with these options.'
691!           stop
692!         endif
693!       endif
694!
695!       if (.not.(doubleq).and.(h2o).and.
696!     $     (chem)) then
697!         if (dustbin.gt.0) then
698!           print*,dustbin,' dust bins'
699!         endif
700!         print*,nqmx-2-dustbin,' chemistry tracers'
701!         print*,' 1 water vapour tracer'
702!         print*,' 1 water ice tracer'
703!       endif
704!
705!       if (.not.(doubleq).and.(h2o).and.
706!     $     .not.(chem)) then
707!         if (dustbin.gt.0) then
708!           print*,dustbin,' dust bins'
709!         endif
710!         print*,' 1 water vapour tracer'
711!         print*,' 1 water ice tracer'
712!         if (nqmx.gt.(dustbin+2)) then
713!           print*,'nqmx should be ',(dustbin+2),
714!     $            ' with these options...'
715!                  print*,'(or check callphys.def)'
716!         endif
717!         if (nqmx.lt.(dustbin+2)) then
718!           write(*,*) "inifis: nqmx.lt.(dustbin+2)"
719!           stop
720!         endif
721!       endif
722!
723!      endif ! of if (tracer)
724!
725!      RETURN
726      END
Note: See TracBrowser for help on using the repository browser.