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

Last change on this file since 278 was 234, checked in by aslmd, 14 years ago

LMDZ.MARS
+ MESOSCALE

  • An important change to merge physiq.F and inifis.F for GCM and mesoscale
  • This is mostly transparent to GCM users and developers (use of MESOSCALE precompiler flags)
  • Makes it easy (and mostly automatic!) for changes in GCM physics to be impacted in mesoscale physics
  • A few minor changes have followed in the GCM (slope scheme, one-point diagnostic).
  • Compilation + run is OK on both sides (GCM and mesoscale).
  • On the mesoscale side, call_meso_physiq?.inc and call_meso_inifis?.inc have been changed accordingly.

Here is the excerpt from README file:

19/07/2011 == AS

  • Finished converging meso_physiq.F and meso_inifis.F towards physiq.F and inifis.F --> see previous point 15/07/2011 --> meso_ routines no longer exist (everything is in meso_inc and transparent to GCM users) --> GCM routines include mesoscale parts within MESOSCALE precompiler commands --> MESOSCALE parts are as hidden as possible not to mess up with GCM users/developers
  • Cleaned inelegant or useless #ifdef [or] #ifndef MESOSCALE in physiq and inifis so that a minimum amount of such precompiler commands is now reached [mainly related to I/O]
  • Added the SF08 slope insolation model in the general physics parameterizations. Added a callslope keyword in inifis.F and callkeys.h --> This keyword is False by default / True if you use -DMESOSCALE
  • Removed the obsolete call to Viking Lander 1 diagnostic Replaced it with a diagnostic for opacity at the domain center [valid for GCM and mesoscale]
File size: 21.6 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         write(*,*) "photochemistry: include chemical species"
416         photochem=.false. ! default value
417         call getin("photochem",photochem)
418         write(*,*) " photochem = ",photochem
419
420
421! THERMOSPHERE
422
423         write(*,*) "call thermosphere ?"
424         callthermos=.false. ! default value
425         call getin("callthermos",callthermos)
426         write(*,*) " callthermos = ",callthermos
427         
428
429         write(*,*) " water included without cycle ",
430     &              "(only if water=.false.)"
431         thermoswater=.false. ! default value
432         call getin("thermoswater",thermoswater)
433         write(*,*) " thermoswater = ",thermoswater
434
435         write(*,*) "call thermal conduction ?",
436     &    " (only if callthermos=.true.)"
437         callconduct=.false. ! default value
438         call getin("callconduct",callconduct)
439         write(*,*) " callconduct = ",callconduct
440
441         write(*,*) "call EUV heating ?",
442     &   " (only if callthermos=.true.)"
443         calleuv=.false.  ! default value
444         call getin("calleuv",calleuv)
445         write(*,*) " calleuv = ",calleuv
446
447         write(*,*) "call molecular viscosity ?",
448     &   " (only if callthermos=.true.)"
449         callmolvis=.false. ! default value
450         call getin("callmolvis",callmolvis)
451         write(*,*) " callmolvis = ",callmolvis
452
453         write(*,*) "call molecular diffusion ?",
454     &   " (only if callthermos=.true.)"
455         callmoldiff=.false. ! default value
456         call getin("callmoldiff",callmoldiff)
457         write(*,*) " callmoldiff = ",callmoldiff
458         
459
460         write(*,*) "call thermospheric photochemistry ?",
461     &   " (only if callthermos=.true.)"
462         thermochem=.false. ! default value
463         call getin("thermochem",thermochem)
464         write(*,*) " thermochem = ",thermochem
465
466         write(*,*) "date for solar flux calculation:",
467     &   " (1985 < date < 2002)"
468         write(*,*) "(Solar min=1996.4 ave=1993.4 max=1990.6)"
469         solarcondate=1993.4 ! default value
470         call getin("solarcondate",solarcondate)
471         write(*,*) " solarcondate = ",solarcondate
472         
473
474         if (.not.callthermos) then
475           if (thermoswater) then
476             print*,'if thermoswater is set, callthermos must be true'
477             stop
478           endif         
479           if (callconduct) then
480             print*,'if callconduct is set, callthermos must be true'
481             stop
482           endif       
483           if (calleuv) then
484             print*,'if calleuv is set, callthermos must be true'
485             stop
486           endif         
487           if (callmolvis) then
488             print*,'if callmolvis is set, callthermos must be true'
489             stop
490           endif       
491           if (callmoldiff) then
492             print*,'if callmoldiff is set, callthermos must be true'
493             stop
494           endif         
495           if (thermochem) then
496             print*,'if thermochem is set, callthermos must be true'
497             stop
498           endif         
499        endif
500
501! Test of incompatibility:
502! if photochem is used, then water should be used too
503
504         if (photochem.and..not.water) then
505           print*,'if photochem is used, water should be used too'
506           stop
507         endif
508
509! if callthermos is used, then thermoswater should be used too
510! (if water not used already)
511
512         if (callthermos .and. .not.water) then
513           if (callthermos .and. .not.thermoswater) then
514             print*,'if callthermos is used, water or thermoswater
515     &               should be used too'
516             stop
517           endif
518         endif
519
520         PRINT*,'--------------------------------------------'
521         PRINT*
522         PRINT*
523      ELSE
524         write(*,*)
525         write(*,*) 'Cannot read file callphys.def. Is it here ?'
526         stop
527      ENDIF
528
5298000  FORMAT(t5,a12,l8)
5308001  FORMAT(t5,a12,i8)
531
532      PRINT*
533      PRINT*,'inifis: daysec',daysec
534      PRINT*
535      PRINT*,'inifis: The radiative transfer is computed:'
536      PRINT*,'           each ',iradia,' physical time-step'
537      PRINT*,'        or each ',iradia*dtphys,' seconds'
538      PRINT*
539! --------------------------------------------------------------
540!  Managing the Longwave radiative transfer
541! --------------------------------------------------------------
542
543!     In most cases, the run just use the following values :
544!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
545      callemis=.true.     
546!     ilwd=10*int(daysec/dtphys) ! bug before 22/10/01       
547      ilwd=1
548      ilwn=1 !2
549      ilwb=1 !2
550      linear=.true.       
551      ncouche=3
552      alphan=0.4
553      semi=0
554
555!     BUT people working hard on the LW may want to read them in 'radia.def'
556!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
557
558      OPEN(99,file='radia.def',status='old',form='formatted'
559     .     ,iostat=ierr)
560      IF(ierr.EQ.0) THEN
561         write(*,*) 'inifis: Reading radia.def !!!'
562         READ(99,fmt='(a)') ch1
563         READ(99,*) callemis
564         WRITE(*,8000) ch1,callemis
565
566         READ(99,fmt='(a)') ch1
567         READ(99,*) iradia
568         WRITE(*,8001) ch1,iradia
569
570         READ(99,fmt='(a)') ch1
571         READ(99,*) ilwd
572         WRITE(*,8001) ch1,ilwd
573
574         READ(99,fmt='(a)') ch1
575         READ(99,*) ilwn
576         WRITE(*,8001) ch1,ilwn
577
578         READ(99,fmt='(a)') ch1
579         READ(99,*) linear
580         WRITE(*,8000) ch1,linear
581
582         READ(99,fmt='(a)') ch1
583         READ(99,*) ncouche
584         WRITE(*,8001) ch1,ncouche
585
586         READ(99,fmt='(a)') ch1
587         READ(99,*) alphan
588         WRITE(*,*) ch1,alphan
589
590         READ(99,fmt='(a)') ch1
591         READ(99,*) ilwb
592         WRITE(*,8001) ch1,ilwb
593
594
595         READ(99,fmt='(a)') ch1
596         READ(99,'(l1)') callg2d
597         WRITE(*,8000) ch1,callg2d
598
599         READ(99,fmt='(a)') ch1
600         READ(99,*) semi
601         WRITE(*,*) ch1,semi
602      end if
603      CLOSE(99)
604
605!-----------------------------------------------------------------------
606!     Some more initialization:
607!     ------------------------
608
609      ! in 'comgeomfi.h'
610      CALL SCOPY(ngrid,plon,1,long,1)
611      CALL SCOPY(ngrid,plat,1,lati,1)
612      CALL SCOPY(ngrid,parea,1,area,1)
613      totarea=SSUM(ngridmx,area,1)
614
615      ! in 'comdiurn.h'
616      DO ig=1,ngrid
617         sinlat(ig)=sin(plat(ig))
618         coslat(ig)=cos(plat(ig))
619         sinlon(ig)=sin(plon(ig))
620         coslon(ig)=cos(plon(ig))
621      ENDDO
622
623      pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h
624
625!     managing the tracers, and tests:
626!     -------------------------------
627!     Ehouarn: removed; as these tests are now done in initracer.F
628!      if(tracer) then
629!
630!!          when photochem is used, nqchem_min is the rank
631!!          of the first chemical species
632!
633!! Ehouarn: nqchem_min is now meaningless and no longer used
634!!       nqchem_min = 1
635!       if (photochem .or. callthermos) then
636!         chem = .true.
637!       end if
638!
639!       if (water .or. thermoswater) h2o = .true.
640!
641!!          TESTS
642!
643!       print*,'inifis: TRACERS:'
644!       write(*,*) "    chem=",chem,"    h2o=",h2o
645!!       write(*,*) "   doubleq=",doubleq
646!!       write(*,*) "   dustbin=",dustbin
647!
648!       if ((doubleq).and.(h2o).and.
649!     $     (chem)) then
650!         print*,' 2 dust tracers (doubleq)'
651!         print*,' 1 water vapour tracer'
652!         print*,' 1 water ice tracer'
653!         print*,nqmx-4,' chemistry tracers'
654!       endif
655!
656!       if ((doubleq).and.(h2o).and.
657!     $     .not.(chem)) then
658!         print*,' 2 dust tracers (doubleq)'
659!         print*,' 1 water vapour tracer'
660!         print*,' 1 water ice tracer'
661!         if (nqmx.LT.4) then
662!           print*,'nqmx should be at least equal to'
663!           print*,'4 with these options.'
664!           stop
665!         endif
666!       endif
667!
668!       if (.not.(doubleq).and.(h2o).and.
669!     $     (chem)) then
670!         if (dustbin.gt.0) then
671!           print*,dustbin,' dust bins'
672!         endif
673!         print*,nqmx-2-dustbin,' chemistry tracers'
674!         print*,' 1 water vapour tracer'
675!         print*,' 1 water ice tracer'
676!       endif
677!
678!       if (.not.(doubleq).and.(h2o).and.
679!     $     .not.(chem)) then
680!         if (dustbin.gt.0) then
681!           print*,dustbin,' dust bins'
682!         endif
683!         print*,' 1 water vapour tracer'
684!         print*,' 1 water ice tracer'
685!         if (nqmx.gt.(dustbin+2)) then
686!           print*,'nqmx should be ',(dustbin+2),
687!     $            ' with these options...'
688!                  print*,'(or check callphys.def)'
689!         endif
690!         if (nqmx.lt.(dustbin+2)) then
691!           write(*,*) "inifis: nqmx.lt.(dustbin+2)"
692!           stop
693!         endif
694!       endif
695!
696!      endif ! of if (tracer)
697!
698!      RETURN
699      END
Note: See TracBrowser for help on using the repository browser.