source: LMDZ.3.3/trunk/libf/phylmd/convect1.F @ 1117

Last change on this file since 1117 was 254, checked in by lmdz, 23 years ago

Inclusion de la version vectorisee de KE FH
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 21.0 KB
Line 
1      subroutine convect1(len,nd,ndp1,noff,minorig,
2     &                   t,q,qs,u,v,
3     &                   p,ph,iflag,ft,
4     &                   fq,fu,fv,precip,cbmf,delt,Ma)
5C.............................START PROLOGUE............................
6C
7C SCCS IDENTIFICATION:  @(#)convect1.f  1.1 04/21/00
8C                       19:40:52 /h/cm/library/nogaps4/src/sub/fcst/convect1.f_v
9C
10C CONFIGURATION IDENTIFICATION:  None
11C
12C MODULE NAME:  convect1
13C
14C DESCRIPTION:
15C
16C convect1     The Emanuel Cumulus Convection Scheme
17C
18C CONTRACT NUMBER AND TITLE:  None
19C
20C REFERENCES: Programmers  K. Emanuel (MIT), Timothy F. Hogan, M. Peng (NRL)
21C
22C CLASSIFICATION:  Unclassified
23C
24C RESTRICTIONS: None
25C
26C COMPILER DEPENDENCIES: FORTRAN 77, FORTRAN 90
27C
28C COMPILE OPTIONS: Fortran 77: -Zu -Wf"-ei -o aggress"
29C                  Fortran 90: -O vector3,scalar3,task1,aggress,overindex  -ei -r 2
30C
31C LIBRARIES OF RESIDENCE: /a/ops/lib/libfcst159.a
32C
33C USAGE: call convect1(len,nd,noff,minorig,
34C    &                   t,q,qs,u,v,
35C    &                   p,ph,iflag,ft,
36C    &                   fq,fu,fv,precip,cbmf,delt)
37C
38C PARAMETERS:
39C      Name            Type         Usage            Description
40C   ----------      ----------     -------  ----------------------------
41C
42C      len           Integer        Input        first (i) dimension
43C      nd            Integer        Input        vertical (k) dimension
44C      ndp1          Integer        Input        nd + 1
45C      noff          Integer        Input        integer limit for convection (nd-noff)
46C      minorig       Integer        Input        First level of convection
47C      t             Real           Input        temperature
48C      q             Real           Input        specific hum
49C      qs            Real           Input        sat specific hum
50C      u             Real           Input        u-wind
51C      v             Real           Input        v-wind
52C      p             Real           Input        full level pressure
53C      ph            Real           Input        half level pressure
54C      iflag         Integer        Output       iflag on latitude strip
55C      ft            Real           Output       temp tend
56C      fq            Real           Output       spec hum tend
57C      fu            Real           Output       u-wind tend
58C      fv            Real           Output       v-wind tend
59C      cbmf          Real           In/Out       cumulus mass flux
60C      delt          Real           Input        time step
61C      iflag         Integer        Output       integer flag for Emanuel conditions
62C
63C COMMON BLOCKS:
64C      Block      Name     Type    Usage              Notes
65C     --------  --------   ----    ------   ------------------------
66C
67C FILES: None
68C
69C DATA BASES: None
70C
71C NON-FILE INPUT/OUTPUT: None
72C
73C ERROR CONDITIONS: None
74C
75C ADDITIONAL COMMENTS: None
76C
77C.................MAINTENANCE SECTION................................
78C
79C MODULES CALLED:
80C         Name           Description
81C         convect2        Emanuel cumulus convection tendency calculations
82C        -------     ----------------------
83C LOCAL VARIABLES AND
84C          STRUCTURES:
85C Name     Type    Description
86C -------  ------  -----------
87C See Comments Below
88C
89C i        Integer loop index
90C k        Integer loop index
91c
92C METHOD:
93C
94C See Emanuel, K. and M. Zivkovic-Rothman, 2000: Development and evaluation of a
95C       convective scheme for use in climate models.
96C
97C FILES: None
98C
99C INCLUDE FILES: None
100C
101C MAKEFILE: /a/ops/met/nogaps/src/sub/fcst/fcst159lib.mak
102C
103C..............................END PROLOGUE.............................
104c
105c
106      implicit none
107c
108#include "dimensions.h"
109#include "dimphy.h"
110c
111      integer len
112      integer nd
113      integer ndp1
114      integer noff
115      real t(len,nd)
116      real q(len,nd)
117      real qs(len,nd)
118      real u(len,nd)
119      real v(len,nd)
120      real p(len,nd)
121      real ph(len,ndp1)
122      integer iflag(len)
123      real ft(len,nd)
124      real fq(len,nd)
125      real fu(len,nd)
126      real fv(len,nd)
127      real precip(len)
128      real cbmf(len)
129      real Ma(len,nd)
130      integer minorig
131      real delt,cpd,cpv,cl,rv,rd,lv0,g
132      real sigs,sigd,elcrit,tlcrit,omtsnow,dtmax,damp
133      real alpha,entp,coeffs,coeffr,omtrain,cu
134c
135!-------------------------------------------------------------------
136! --- ARGUMENTS
137!-------------------------------------------------------------------
138! --- On input:
139!
140!  t:   Array of absolute temperature (K) of dimension ND, with first
141!       index corresponding to lowest model level. Note that this array
142!       will be altered by the subroutine if dry convective adjustment
143!       occurs and if IPBL is not equal to 0.
144!
145!  q:   Array of specific humidity (gm/gm) of dimension ND, with first
146!       index corresponding to lowest model level. Must be defined
147!       at same grid levels as T. Note that this array will be altered
148!       if dry convective adjustment occurs and if IPBL is not equal to 0.
149!
150!  qs:  Array of saturation specific humidity of dimension ND, with first
151!       index corresponding to lowest model level. Must be defined
152!       at same grid levels as T. Note that this array will be altered
153!       if dry convective adjustment occurs and if IPBL is not equal to 0.
154!
155!  u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
156!       index corresponding with the lowest model level. Defined at
157!       same levels as T. Note that this array will be altered if
158!       dry convective adjustment occurs and if IPBL is not equal to 0.
159!
160!  v:   Same as u but for meridional velocity.
161!
162!  tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
163!       where NTRA is the number of different tracers. If no
164!       convective tracer transport is needed, define a dummy
165!       input array of dimension (ND,1). Tracers are defined at
166!       same vertical levels as T. Note that this array will be altered
167!       if dry convective adjustment occurs and if IPBL is not equal to 0.
168!
169!  p:   Array of pressure (mb) of dimension ND, with first
170!       index corresponding to lowest model level. Must be defined
171!       at same grid levels as T.
172!
173!  ph:  Array of pressure (mb) of dimension ND+1, with first index
174!       corresponding to lowest level. These pressures are defined at
175!       levels intermediate between those of P, T, Q and QS. The first
176!       value of PH should be greater than (i.e. at a lower level than)
177!       the first value of the array P.
178!
179!  nl:  The maximum number of levels to which convection can penetrate, plus 1.
180!       NL MUST be less than or equal to ND-1.
181!
182!  delt: The model time step (sec) between calls to CONVECT
183!
184!----------------------------------------------------------------------------
185! ---   On Output:
186!
187!  iflag: An output integer whose value denotes the following:
188!       VALUE   INTERPRETATION
189!       -----   --------------
190!         0     Moist convection occurs.
191!         1     Moist convection occurs, but a CFL condition
192!               on the subsidence warming is violated. This
193!               does not cause the scheme to terminate.
194!         2     Moist convection, but no precip because ep(inb) lt 0.0001
195!         3     No moist convection because new cbmf is 0 and old cbmf is 0.
196!         4     No moist convection; atmosphere is not
197!               unstable
198!         6     No moist convection because ihmin le minorig.
199!         7     No moist convection because unreasonable
200!               parcel level temperature or specific humidity.
201!         8     No moist convection: lifted condensation
202!               level is above the 200 mb level.
203!         9     No moist convection: cloud base is higher
204!               then the level NL-1.
205!
206!  ft:   Array of temperature tendency (K/s) of dimension ND, defined at same
207!        grid levels as T, Q, QS and P.
208!
209!  fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
210!        defined at same grid levels as T, Q, QS and P.
211!
212!  fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
213!        defined at same grid levels as T.
214!
215!  fv:   Same as FU, but for forcing of meridional velocity.
216!
217!  ftra: Array of forcing of tracer content, in tracer mixing ratio per
218!        second, defined at same levels as T. Dimensioned (ND,NTRA).
219!
220!  precip: Scalar convective precipitation rate (mm/day).
221!
222!  wd:   A convective downdraft velocity scale. For use in surface
223!        flux parameterizations. See convect.ps file for details.
224!
225!  tprime: A convective downdraft temperature perturbation scale (K).
226!          For use in surface flux parameterizations. See convect.ps
227!          file for details.
228!
229!  qprime: A convective downdraft specific humidity
230!          perturbation scale (gm/gm).
231!          For use in surface flux parameterizations. See convect.ps
232!          file for details.
233!
234!  cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
235!        BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
236!        ITS NEXT CALL. That is, the value of CBMF must be "remembered"
237!        by the calling program between calls to CONVECT.
238!
239!  det:   Array of detrainment mass flux of dimension ND.
240!
241!-------------------------------------------------------------------
242c
243c  Local arrays
244c
245      integer nl
246      integer nlp
247      integer nlm
248      integer i,k,n
249      real delti
250      real rowl
251      real clmcpv
252      real clmcpd
253      real cpdmcp
254      real cpvmcpd
255      real eps
256      real epsi
257      real epsim1
258      real ginv
259      real hrd
260      real prccon1
261      integer icbmax
262      real lv(klon,klev)
263      real cpn(klon,klev)
264      real cpx(klon,klev)
265      real tv(klon,klev)
266      real gz(klon,klev)
267      real hm(klon,klev)
268      real h(klon,klev)
269      real work(klon)
270      integer ihmin(klon)
271      integer nk(klon)
272      real rh(klon)
273      real chi(klon)
274      real plcl(klon)
275      integer icb(klon)
276      real tnk(klon)
277      real qnk(klon)
278      real gznk(klon)
279      real pnk(klon)
280      real qsnk(klon)
281      real ticb(klon)
282      real gzicb(klon)
283      real tp(klon,klev)
284      real tvp(klon,klev)
285      real clw(klon,klev)
286c
287      real ah0(klon),cpp(klon)
288      real tg,qg,s,alv,tc,ahg,denom,es,rg
289c
290      integer ncum
291      integer idcum(klon)
292c
293      cpd=1005.7
294      cpv=1870.0
295      cl=4190.0
296      rv=461.5
297      rd=287.04
298      lv0=2.501E6
299      g=9.8
300C
301C   *** ELCRIT IS THE AUTOCONVERSION THERSHOLD WATER CONTENT (gm/gm) ***
302C   ***  TLCRIT IS CRITICAL TEMPERATURE BELOW WHICH THE AUTO-        ***
303C   ***       CONVERSION THRESHOLD IS ASSUMED TO BE ZERO             ***
304C   ***     (THE AUTOCONVERSION THRESHOLD VARIES LINEARLY            ***
305C   ***               BETWEEN 0 C AND TLCRIT)                        ***
306C   ***   ENTP IS THE COEFFICIENT OF MIXING IN THE ENTRAINMENT       ***
307C   ***                       FORMULATION                            ***
308C   ***  SIGD IS THE FRACTIONAL AREA COVERED BY UNSATURATED DNDRAFT  ***
309C   ***  SIGS IS THE FRACTION OF PRECIPITATION FALLING OUTSIDE       ***
310C   ***                        OF CLOUD                              ***
311C   ***        OMTRAIN IS THE ASSUMED FALL SPEED (P/s) OF RAIN       ***
312C   ***     OMTSNOW IS THE ASSUMED FALL SPEED (P/s) OF SNOW          ***
313C   ***  COEFFR IS A COEFFICIENT GOVERNING THE RATE OF EVAPORATION   ***
314C   ***                          OF RAIN                             ***
315C   ***  COEFFS IS A COEFFICIENT GOVERNING THE RATE OF EVAPORATION   ***
316C   ***                          OF SNOW                             ***
317C   ***     CU IS THE COEFFICIENT GOVERNING CONVECTIVE MOMENTUM      ***
318C   ***                         TRANSPORT                            ***
319C   ***    DTMAX IS THE MAXIMUM NEGATIVE TEMPERATURE PERTURBATION    ***
320C   ***        A LIFTED PARCEL IS ALLOWED TO HAVE BELOW ITS LFC      ***
321C   ***    ALPHA AND DAMP ARE PARAMETERS THAT CONTROL THE RATE OF    ***
322C   ***                 APPROACH TO QUASI-EQUILIBRIUM                ***
323C   ***   (THEIR STANDARD VALUES ARE  0.20 AND 0.1, RESPECTIVELY)    ***
324C   ***                   (DAMP MUST BE LESS THAN 1)                 ***
325c
326      sigs=0.12
327      sigd=0.05
328      elcrit=0.0011
329      tlcrit=-55.0
330      omtsnow=5.5
331      dtmax=0.9
332      damp=0.1
333      alpha=0.2
334      entp=1.5
335      coeffs=0.8
336      coeffr=1.0
337      omtrain=50.0
338c
339      cu=0.70
340      damp=0.1
341c
342c
343c Define nl, nlp, nlm, and delti
344c
345      nl=nd-noff
346      nlp=nl+1
347      nlm=nl-1
348      delti=1.0/delt
349!
350!-------------------------------------------------------------------
351! --- SET CONSTANTS
352!-------------------------------------------------------------------
353!
354      rowl=1000.0
355      clmcpv=cl-cpv
356      clmcpd=cl-cpd
357      cpdmcp=cpd-cpv
358      cpvmcpd=cpv-cpd
359      eps=rd/rv
360      epsi=1.0/eps
361      epsim1=epsi-1.0
362      ginv=1.0/g
363      hrd=0.5*rd
364      prccon1=86400.0*1000.0/(rowl*g)
365!
366! dtmax is the maximum negative temperature perturbation.
367!
368!=====================================================================
369! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
370!=====================================================================
371!
372      do 20 k=1,nd
373        do 10 i=1,len
374         ft(i,k)=0.0
375         fq(i,k)=0.0
376         fu(i,k)=0.0
377         fv(i,k)=0.0
378         tvp(i,k)=0.0
379         tp(i,k)=0.0
380         clw(i,k)=0.0
381         gz(i,k) = 0.
382 10     continue
383 20   continue
384      do 60 i=1,len
385        precip(i)=0.0
386        iflag(i)=0
387 60   continue
388c
389!=====================================================================
390! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
391!=====================================================================
392      do 110 k=1,nl+1
393        do 100 i=1,len
394          lv(i,k)= lv0-clmcpv*(t(i,k)-273.15)
395          cpn(i,k)=cpd*(1.0-q(i,k))+cpv*q(i,k)
396          cpx(i,k)=cpd*(1.0-q(i,k))+cl*q(i,k)
397          tv(i,k)=t(i,k)*(1.0+q(i,k)*epsim1)
398 100    continue
399 110  continue
400c
401c gz = phi at the full levels (same as p).
402c
403      do 120 i=1,len
404        gz(i,1)=0.0
405 120  continue
406      do 140 k=2,nlp
407        do 130 i=1,len
408          gz(i,k)=gz(i,k-1)+hrd*(tv(i,k-1)+tv(i,k))
409     &         *(p(i,k-1)-p(i,k))/ph(i,k)
410 130    continue
411 140  continue
412c
413c h  = phi + cpT (dry static energy).
414c hm = phi + cp(T-Tbase)+Lq
415c
416      do 170 k=1,nlp
417        do 160 i=1,len
418          h(i,k)=gz(i,k)+cpn(i,k)*t(i,k)
419          hm(i,k)=gz(i,k)+cpx(i,k)*(t(i,k)-t(i,1))+lv(i,k)*q(i,k)
420 160    continue
421 170  continue
422c
423!-------------------------------------------------------------------
424! --- Find level of minimum moist static energy
425! --- If level of minimum moist static energy coincides with
426! --- or is lower than minimum allowable parcel origin level,
427! --- set iflag to 6.
428!-------------------------------------------------------------------
429      do 180 i=1,len
430       work(i)=1.0e12
431       ihmin(i)=nl
432 180  continue
433      do 200 k=2,nlp
434        do 190 i=1,len
435         if((hm(i,k).lt.work(i)).and.
436     &      (hm(i,k).lt.hm(i,k-1)))then
437           work(i)=hm(i,k)
438           ihmin(i)=k
439         endif
440 190    continue
441 200  continue
442      do 210 i=1,len
443        ihmin(i)=min(ihmin(i),nlm)
444        if(ihmin(i).le.minorig)then
445          iflag(i)=6
446        endif
447 210  continue
448c
449!-------------------------------------------------------------------
450! --- Find that model level below the level of minimum moist static
451! --- energy that has the maximum value of moist static energy
452!-------------------------------------------------------------------
453 
454      do 220 i=1,len
455       work(i)=hm(i,minorig)
456       nk(i)=minorig
457 220  continue
458      do 240 k=minorig+1,nl
459        do 230 i=1,len
460         if((hm(i,k).gt.work(i)).and.(k.le.ihmin(i)))then
461           work(i)=hm(i,k)
462           nk(i)=k
463         endif
464 230     continue
465 240  continue
466!-------------------------------------------------------------------
467! --- Check whether parcel level temperature and specific humidity
468! --- are reasonable
469!-------------------------------------------------------------------
470       do 250 i=1,len
471       if(((t(i,nk(i)).lt.250.0).or.
472     &      (q(i,nk(i)).le.0.0).or.
473     &      (p(i,ihmin(i)).lt.400.0)).and.
474     &      (iflag(i).eq.0))iflag(i)=7
475 250   continue
476!-------------------------------------------------------------------
477! --- Calculate lifted condensation level of air at parcel origin level
478! --- (Within 0.2% of formula of Bolton, MON. WEA. REV.,1980)
479!-------------------------------------------------------------------
480       do 260 i=1,len
481        tnk(i)=t(i,nk(i))
482        qnk(i)=q(i,nk(i))
483        gznk(i)=gz(i,nk(i))
484        pnk(i)=p(i,nk(i))
485        qsnk(i)=qs(i,nk(i))
486c
487        rh(i)=qnk(i)/qsnk(i)
488        rh(i)=min(1.0,rh(i))
489        chi(i)=tnk(i)/(1669.0-122.0*rh(i)-tnk(i))
490        plcl(i)=pnk(i)*(rh(i)**chi(i))
491        if(((plcl(i).lt.200.0).or.(plcl(i).ge.2000.0))
492     &   .and.(iflag(i).eq.0))iflag(i)=8
493 260   continue
494!-------------------------------------------------------------------
495! --- Calculate first level above lcl (=icb)
496!-------------------------------------------------------------------
497      do 270 i=1,len
498       icb(i)=nlm
499 270  continue
500c
501      do 290 k=minorig,nl
502        do 280 i=1,len
503          if((k.ge.(nk(i)+1)).and.(p(i,k).lt.plcl(i)))
504     &    icb(i)=min(icb(i),k)
505 280    continue
506 290  continue
507c
508      do 300 i=1,len
509        if((icb(i).ge.nlm).and.(iflag(i).eq.0))iflag(i)=9
510 300  continue
511c
512c Compute icbmax.
513c
514      icbmax=2
515      do 310 i=1,len
516        icbmax=max(icbmax,icb(i))
517 310  continue
518!
519!-------------------------------------------------------------------
520! --- Calculates the lifted parcel virtual temperature at nk,
521! --- the actual temperature, and the adiabatic
522! --- liquid water content. The procedure is to solve the equation.
523!     cp*tp+L*qp+phi=cp*tnk+L*qnk+gznk.
524!-------------------------------------------------------------------
525!
526      do 320 i=1,len
527        tnk(i)=t(i,nk(i))
528        qnk(i)=q(i,nk(i))
529        gznk(i)=gz(i,nk(i))
530        ticb(i)=t(i,icb(i))
531        gzicb(i)=gz(i,icb(i))
532 320  continue
533c
534c   ***  Calculate certain parcel quantities, including static energy   ***
535c
536      do 330 i=1,len
537        ah0(i)=(cpd*(1.-qnk(i))+cl*qnk(i))*tnk(i)
538     &         +qnk(i)*(lv0-clmcpv*(tnk(i)-273.15))+gznk(i)
539        cpp(i)=cpd*(1.-qnk(i))+qnk(i)*cpv
540 330  continue
541c
542c   ***   Calculate lifted parcel quantities below cloud base   ***
543c
544        do 350 k=minorig,icbmax-1
545          do 340 i=1,len
546           tp(i,k)=tnk(i)-(gz(i,k)-gznk(i))/cpp(i)
547           tvp(i,k)=tp(i,k)*(1.+qnk(i)*epsi)
548  340     continue
549  350   continue
550c
551c    ***  Find lifted parcel quantities above cloud base    ***
552c
553        do 360 i=1,len
554         tg=ticb(i)
555         qg=qs(i,icb(i))
556         alv=lv0-clmcpv*(ticb(i)-273.15)
557c
558c First iteration.
559c
560          s=cpd+alv*alv*qg/(rv*ticb(i)*ticb(i))
561          s=1./s
562          ahg=cpd*tg+(cl-cpd)*qnk(i)*ticb(i)+alv*qg+gzicb(i)
563          tg=tg+s*(ah0(i)-ahg)
564          tg=max(tg,35.0)
565          tc=tg-273.15
566          denom=243.5+tc
567          if(tc.ge.0.0)then
568           es=6.112*exp(17.67*tc/denom)
569          else
570           es=exp(23.33086-6111.72784/tg+0.15215*log(tg))
571          endif
572          qg=eps*es/(p(i,icb(i))-es*(1.-eps))
573c
574c Second iteration.
575c
576          s=cpd+alv*alv*qg/(rv*ticb(i)*ticb(i))
577          s=1./s
578          ahg=cpd*tg+(cl-cpd)*qnk(i)*ticb(i)+alv*qg+gzicb(i)
579          tg=tg+s*(ah0(i)-ahg)
580          tg=max(tg,35.0)
581          tc=tg-273.15
582          denom=243.5+tc
583          if(tc.ge.0.0)then
584           es=6.112*exp(17.67*tc/denom)
585          else
586           es=exp(23.33086-6111.72784/tg+0.15215*log(tg))
587          end if
588          qg=eps*es/(p(i,icb(i))-es*(1.-eps))
589c
590         alv=lv0-clmcpv*(ticb(i)-273.15)
591         tp(i,icb(i))=(ah0(i)-(cl-cpd)*qnk(i)*ticb(i)
592     &   -gz(i,icb(i))-alv*qg)/cpd
593         clw(i,icb(i))=qnk(i)-qg
594         clw(i,icb(i))=max(0.0,clw(i,icb(i)))
595         rg=qg/(1.-qnk(i))
596         tvp(i,icb(i))=tp(i,icb(i))*(1.+rg*epsi)
597  360   continue
598c
599      do 380 k=minorig,icbmax
600       do 370 i=1,len
601         tvp(i,k)=tvp(i,k)-tp(i,k)*qnk(i)
602 370   continue
603 380  continue
604c
605!-------------------------------------------------------------------
606! --- Test for instability.
607! --- If there was no convection at last time step and parcel
608! --- is stable at icb, then set iflag to 4.
609!-------------------------------------------------------------------
610 
611      do 390 i=1,len
612        if((cbmf(i).eq.0.0) .and.(iflag(i).eq.0).and.
613     &  (tvp(i,icb(i)).le.(tv(i,icb(i))-dtmax)))iflag(i)=4
614 390  continue
615 
616!=====================================================================
617! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
618!=====================================================================
619c
620      ncum=0
621      do 400 i=1,len
622        if(iflag(i).eq.0)then
623           ncum=ncum+1
624           idcum(ncum)=i
625        endif
626 400  continue
627c
628c Call convect2, which compresses the points and computes the heating,
629c moistening, velocity mixing, and precipiation.
630c
631c     print*,'cpd avant convect2 ',cpd
632      if(ncum.gt.0)then
633      call convect2(ncum,idcum,len,nd,ndp1,nl,minorig,
634     &              nk,icb,
635     &              t,q,qs,u,v,gz,tv,tp,tvp,clw,h,
636     &              lv,cpn,p,ph,ft,fq,fu,fv,
637     &              tnk,qnk,gznk,plcl,
638     &              precip,cbmf,iflag,
639     &              delt,cpd,cpv,cl,rv,rd,lv0,g,
640     &              sigs,sigd,elcrit,tlcrit,omtsnow,dtmax,damp,
641     &              alpha,entp,coeffs,coeffr,omtrain,cu,Ma)
642      endif
643c
644      return
645      end
Note: See TracBrowser for help on using the repository browser.