56       subroutine sascnvn(im,ix,km,jcap,delt,delp,prslp,psp,phil,ql,     &
 
   57      &     q1,t1,u1,v1,cldwrk,rn,kbot,ktop,kcnv,islimsk,                &
 
   58      &     dot,ncloud,ud_mf,dd_mf,dt_mf,cnvw,cnvc)
 
   62       use machine 
, only : kind_phys
 
   63       use funcphys 
, only : fpvs
 
   67      &,             eps => con_eps, epsm1 => con_epsm1
 
   70       integer            im, ix,  km, jcap, ncloud,                     &
 
   71      &                   kbot(im), ktop(im), kcnv(im)                  
 
   73       real(kind=kind_phys) delt
 
   74       real(kind=kind_phys) psp(im),    delp(ix,km), prslp(ix,km)
 
   75       real(kind=kind_phys) ps(im),     del(ix,km),  prsl(ix,km),        &
 
   76      &                     ql(ix,km,2),q1(ix,km),   t1(ix,km),          &
 
   77      &                     u1(ix,km),  v1(ix,km),                       & 
 
   78      &                     cldwrk(im), rn(im),                          &
 
   79      &                     dot(ix,km), phil(ix,km),                     &
 
   80      &                     cnvw(ix,km), cnvc(ix,km),                    &
 
   81      &                     ud_mf(im,km),dd_mf(im,km),dt_mf(im,km) 
 
   83       integer              i, indx, jmn, k, kk, km1
 
   84       integer, 
dimension(im), 
intent(in) :: islimsk
 
   87       real(kind=kind_phys) clam, cxlamu, xlamde, xlamdd
 
   90       real(kind=kind_phys) adw,     aup,     aafac,
 
   95      &                     dq,      dqsdp,   dqsdt,   dt,
 
   96      &                     dt2,     dtmax,   dtmin,   dv1h,
 
   97      &                     dv1q,    dv2h,    dv2q,    dv1u,
 
   98      &                     dv1v,    dv2u,    dv2v,    dv3q,
 
  100      &                     dz,      dz1,     e1,      edtmax,
 
  101      &                     edtmaxl, edtmaxs, el2orc,  elocp,
 
  102      &                     es,      etah,    cthk,    dthk,
 
  103      &                     evef,    evfact,  evfactl, fact1,
 
  104      &                     fact2,   factor,  fjcap,   fkm,
 
  107      &                     rain,    rfact,   shear,   tem1,
 
  109      &                     val2,    w1,      w1l,     w1s,
 
  112      &                     w4s,     xdby,    xpw,     xpwd,
 
  116       integer              kb(im), kbcon(im), kbcon1(im),
 
  117      &                     ktcon(im), ktcon1(im),
 
  118      &                     jmin(im), lmin(im), kbmax(im),
 
  121       real(kind=kind_phys) aa1(im),     acrt(im),   acrtfct(im),
 
  122      &                     delhbar(im), delq(im),   delq2(im),
 
  123      &                     delqbar(im), delqev(im), deltbar(im),
 
  124      &                     deltv(im),   dtconv(im), edt(im),
 
  125      &                     edto(im),    edtx(im),   fld(im),
 
  126      &                     hcdo(im,km), hmax(im),   hmin(im),
 
  127      &                     ucdo(im,km), vcdo(im,km),aa2(im),
 
  128      &                     pbcdif(im),  pdot(im),   po(im,km),
 
  129      &                     pwavo(im),   pwevo(im),  xlamud(im),
 
  130      &                     qcdo(im,km), qcond(im),  qevap(im),
 
  131      &                     rntot(im),   vshear(im), xaa0(im),
 
  133      &                     xmb(im),     xmbmax(im), xpwav(im),
 
  134      &                     xpwev(im),   delubar(im),delvbar(im)
 
  136       real(kind=kind_phys) cincr, cincrmax, cincrmin
 
  138 c  physical parameters
 
  140       parameter(elocp=hvap/cp,
 
  141      &          el2orc=hvap*hvap/(rv*cp))
 
  142       parameter(c0=.002,c1=.002,delta=fv)
 
  143       parameter(fact1=(cvap-cliq)/rv,fact2=hvap/rv-fact1*t0c)
 
  144       parameter(cthk=150.,cincrmax=180.,cincrmin=120.,dthk=25.)
 
  145 c  local variables and arrays
 
  146       real(kind=kind_phys) pfld(im,km),    to(im,km),     qo(im,km),
 
  147      &                     uo(im,km),      vo(im,km),     qeso(im,km)
 
  150       real(kind=kind_phys) qlko_ktcon(im), dellal(im,km), tvo(im,km),
 
  151      &                     dbyo(im,km),    zo(im,km),     xlamue(im,km),
 
  152      &                     fent1(im,km),   fent2(im,km),  frh(im,km),
 
  153      &                     heo(im,km),     heso(im,km),
 
  154      &                     qrcd(im,km),    dellah(im,km), dellaq(im,km),
 
  155      &                     dellau(im,km),  dellav(im,km), hcko(im,km),
 
  156      &                     ucko(im,km),    vcko(im,km),   qcko(im,km),
 
  157      &                     eta(im,km),     etad(im,km),   zi(im,km),
 
  158      &                     qrcko(im,km),   qrcdo(im,km),
 
  159      &                     pwo(im,km),     pwdo(im,km),
 
  160      &                     tx1(im),        sumx(im),      cnvwt(im,km)
 
  163       logical totflg, cnvflg(im), flg(im)
 
  165       real(kind=kind_phys) pcrit(15), acritt(15), acrit(15)
 
  167       data pcrit/850.,800.,750.,700.,650.,600.,550.,500.,450.,400.,
 
  168      &           350.,300.,250.,200.,150./
 
  169       data acritt/.0633,.0445,.0553,.0664,.075,.1082,.1521,.2216,
 
  170      &           .3151,.3677,.41,.5255,.7663,1.1686,1.6851/
 
  172 c     
data acritt/.203,.515,.521,.566,.625,.665,.659,.688,
 
  173 c    &            .743,.813,.886,.947,1.138,1.377,1.896/
 
  174       real(kind=kind_phys) tf, tcr, tcrf
 
  175       parameter(tf=233.16, tcr=263.16, tcrf=1.0/(tcr-tf))
 
  177 c-----------------------------------------------------------------------
 
  240         acrit(k) = acritt(k) * (975. - pcrit(k))
 
  244       dtmin = max(dt2, val )
 
  246       dtmax = max(dt2, val )
 
  247 c  model tunable parameters are all here
 
  267       fjcap   = (float(jcap) / 126.) ** 2
 
  269       fjcap   = max(fjcap,val)
 
  270       fkm     = (float(km) / 28.) ** 2
 
  282 c  define top layer for search of the downdraft originating layer
 
  283 c  and the maximum thetae for updraft
 
  294           if (prsl(i,k)*tx1(i) .gt. 0.04) kmax(i)  = k + 1
 
  295           if (prsl(i,k)*tx1(i) .gt. 0.45) kbmax(i) = k + 1
 
  296           if (prsl(i,k)*tx1(i) .gt. 0.70) kbm(i)   = k + 1
 
  300         kmax(i)  = min(km,kmax(i))
 
  301         kbmax(i) = min(kbmax(i),kmax(i))
 
  302         kbm(i)   = min(kbm(i),kmax(i))
 
  305 c  hydrostatic height assume zero terr and initially assume
 
  306 c    updraft entrainment rate as an inverse
 function of height
 
  311           zo(i,k) = phil(i,k) / g
 
  317           zi(i,k) = 0.5*(zo(i,k)+zo(i,k+1))
 
  318           xlamue(i,k) = clam / zi(i,k)
 
  324 c   convert surface pressure to mb from cb
 
  328           if (k .le. kmax(i)) 
then 
  329             pfld(i,k) = prsl(i,k) * 10.0
 
  362 c  p is pressure of the layer(mb)
 
  363 c  t is temperature at t-dt(k)..tn
 
  364 c  q is mixing ratio at t-dt(kg/kg)..qn
 
  365 c  to is temperature at t+dt(k)... this is after advection and turbulan
 
  366 c  qo is mixing ratio at t+dt(kg/kg)..q1
 
  371           if (k .le. kmax(i)) 
then 
  372             qeso(i,k) = 0.01 * fpvs(to(i,k))      
 
  373             qeso(i,k) = eps * qeso(i,k) / (pfld(i,k) + epsm1*qeso(i,k))
 
  375             qeso(i,k) = max(qeso(i,k), val1)
 
  377             qo(i,k)   = max(qo(i,k), val2 )
 
  384 c  compute moist static energy
 
  389           if (k .le. kmax(i)) 
then 
  391             tem       = phil(i,k) + cp * to(i,k)
 
  392             heo(i,k)  = tem  + hvap * qo(i,k)
 
  393             heso(i,k) = tem  + hvap * qeso(i,k)
 
  394 c           heo(i,k)  = min(heo(i,k),heso(i,k))
 
  401 c  determine level with largest moist static energy
 
  402 c  this is the level 
where updraft starts
 
  411           if (k .le. kbm(i)) 
then 
  412             if(heo(i,k).gt.hmax(i)) 
then 
  423           if (k .le. kmax(i)-1) 
then 
  424             dz      = .5 * (zo(i,k+1) - zo(i,k))
 
  425             dp      = .5 * (pfld(i,k+1) - pfld(i,k))
 
  426             es      = 0.01 * fpvs(to(i,k+1))      
 
  427             pprime  = pfld(i,k+1) + epsm1 * es
 
  428             qs      = eps * es / pprime
 
  429             dqsdp   = - qs / pprime
 
  430             desdt   = es * (fact1 / to(i,k+1) + fact2 / (to(i,k+1)**2))
 
  431             dqsdt   = qs * pfld(i,k+1) * desdt / (es * pprime)
 
  432             gamma   = el2orc * qeso(i,k+1) / (to(i,k+1)**2)
 
  433             dt      = (g * dz + hvap * dqsdp * dp) / (cp * (1. + gamma))
 
  434             dq      = dqsdt * dt + dqsdp * dp
 
  435             to(i,k) = to(i,k+1) + dt
 
  436             qo(i,k) = qo(i,k+1) + dq
 
  437             po(i,k) = .5 * (pfld(i,k) + pfld(i,k+1))
 
  445           if (k .le. kmax(i)-1) 
then 
  446             qeso(i,k) = 0.01 * fpvs(to(i,k))      
 
  447             qeso(i,k) = eps * qeso(i,k) / (po(i,k) + epsm1*qeso(i,k))
 
  449             qeso(i,k) = max(qeso(i,k), val1)
 
  451             qo(i,k)   = max(qo(i,k), val2 )
 
  453             frh(i,k)  = 1. - min(qo(i,k)/qeso(i,k), 1.)
 
  454             heo(i,k)  = .5 * g * (zo(i,k) + zo(i,k+1)) +
 
  455      &                  cp * to(i,k) + hvap * qo(i,k)
 
  456             heso(i,k) = .5 * g * (zo(i,k) + zo(i,k+1)) +
 
  457      &                  cp * to(i,k) + hvap * qeso(i,k)
 
  458             uo(i,k)   = .5 * (uo(i,k) + uo(i,k+1))
 
  459             vo(i,k)   = .5 * (vo(i,k) + vo(i,k+1))
 
  464 c  look for the level of free convection as cloud base
 
  473           if (flg(i).and.k.le.kbmax(i)) 
then 
  474             if(k.gt.kb(i).and.heo(i,kb(i)).gt.heso(i,k)) 
then 
  484         if(kbcon(i).eq.kmax(i)) cnvflg(i) = .false.
 
  489         totflg = totflg .and. (.not. cnvflg(i))
 
  494 c  determine critical convective inhibition
 
  495 c  as a
 function of vertical velocity at cloud base.
 
  501           pdot(i)  = 0.01 * dot(i,kbcon(i)) 
 
  506           if(islimsk(i) == 1) 
then 
  517           if(pdot(i).le.w4) 
then 
  518             tem = (pdot(i) - w4) / (w3 - w4)
 
  519           elseif(pdot(i).ge.-w4) 
then 
  520             tem = - (pdot(i) + w4) / (w4 - w3)
 
  529           tem1= .5*(cincrmax-cincrmin)
 
  530           cincr = cincrmax - tem * tem1
 
  531           pbcdif(i) = pfld(i,kb(i)) - pfld(i,kbcon(i))
 
  532           if(pbcdif(i).gt.cincr) 
then 
  540         totflg = totflg .and. (.not. cnvflg(i))
 
  545 c  assume that updraft entrainment rate above cloud base is
 
  546 c    same as that at cloud base
 
  556      &      (k.gt.kbcon(i).and.k.lt.kmax(i))) 
then 
  557               xlamue(i,k) = xlamue(i,kbcon(i))
 
  562 c  assume the detrainment rate for the updrafts to be same as
 
  563 c  the entrainment rate at cloud base
 
  568           xlamud(i) = xlamue(i,kbcon(i))
 
  572 c  functions rapidly decreasing with height, mimicking a cloud ensemble
 
  573 c    (bechtold et al., 2008)
 
  578      &      (k.gt.kbcon(i).and.k.lt.kmax(i))) 
then 
  579               tem = qeso(i,k)/qeso(i,kbcon(i))
 
  586 c  final entrainment rate as the sum of turbulent part and organized entrainment
 
  587 c    depending on the environmental relative humidity
 
  588 c    (bechtold et al., 2008)
 
  593      &      (k.ge.kbcon(i).and.k.lt.kmax(i))) 
then 
  594               tem = cxlamu * frh(i,k) * fent2(i,k)
 
  595               xlamue(i,k) = xlamue(i,k)*fent1(i,k) + tem
 
  600 c  determine updraft mass flux for the subcloud layers
 
  610             if(k.lt.kbcon(i).and.k.ge.kb(i)) 
then 
  611               dz       = zi(i,k+1) - zi(i,k)
 
  612               ptem     = 0.5*(xlamue(i,k)+xlamue(i,k+1))-xlamud(i)
 
  613               eta(i,k) = eta(i,k+1) / (1. + ptem * dz)
 
  619 c  compute mass flux above cloud base
 
  624            if(k.gt.kbcon(i).and.k.lt.kmax(i)) 
then 
  625               dz       = zi(i,k) - zi(i,k-1)
 
  626               ptem     = 0.5*(xlamue(i,k)+xlamue(i,k-1))-xlamud(i)
 
  627               eta(i,k) = eta(i,k-1) * (1 + ptem * dz)
 
  633 c  compute updraft cloud properties
 
  639           hcko(i,indx) = heo(i,indx)
 
  640           ucko(i,indx) = uo(i,indx)
 
  641           vcko(i,indx) = vo(i,indx)
 
  646 c  cloud property is modified by the entrainment process
 
  652             if(k.gt.kb(i).and.k.lt.kmax(i)) 
then 
  653               dz   = zi(i,k) - zi(i,k-1)
 
  654               tem  = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
 
  655               tem1 = 0.5 * xlamud(i) * dz
 
  656               factor = 1. + tem - tem1
 
  657               ptem = 0.5 * tem + pgcon
 
  658               ptem1= 0.5 * tem - pgcon
 
  659               hcko(i,k) = ((1.-tem1)*hcko(i,k-1)+tem*0.5*
 
  660      &                     (heo(i,k)+heo(i,k-1)))/factor
 
  661               ucko(i,k) = ((1.-tem1)*ucko(i,k-1)+ptem*uo(i,k)
 
  662      &                     +ptem1*uo(i,k-1))/factor
 
  663               vcko(i,k) = ((1.-tem1)*vcko(i,k-1)+ptem*vo(i,k)
 
  664      &                     +ptem1*vo(i,k-1))/factor
 
  665               dbyo(i,k) = hcko(i,k) - heso(i,k)
 
  671 c   taking account into convection inhibition due to existence of
 
  672 c    dry layers below cloud base
 
  681         if (flg(i).and.k.lt.kmax(i)) 
then 
  682           if(k.ge.kbcon(i).and.dbyo(i,k).gt.0.) 
then 
  691           if(kbcon1(i).eq.kmax(i)) cnvflg(i) = .false.
 
  696           tem = pfld(i,kbcon(i)) - pfld(i,kbcon1(i))
 
  705         totflg = totflg .and. (.not. cnvflg(i))
 
  710 c  determine first guess cloud top as the level of zero buoyancy
 
  719         if (flg(i).and.k .lt. kmax(i)) 
then 
  720           if(k.gt.kbcon1(i).and.dbyo(i,k).lt.0.) 
then 
  730           tem = pfld(i,kbcon(i))-pfld(i,ktcon(i))
 
  731           if(tem.lt.cthk) cnvflg(i) = .false.
 
  737         totflg = totflg .and. (.not. cnvflg(i))
 
  742 c  search for downdraft originating level above theta-e minimum
 
  747            hmin(i) = heo(i,kbcon1(i))
 
  754           if (cnvflg(i) .and. k .le. kbmax(i)) 
then 
  755             if(k.gt.kbcon1(i).and.heo(i,k).lt.hmin(i)) 
then 
  763 c  make sure that jmin(i) is within the cloud
 
  767           jmin(i) = min(lmin(i),ktcon(i)-1)
 
  768           jmin(i) = max(jmin(i),kbcon1(i)+1)
 
  769           if(jmin(i).ge.ktcon(i)) cnvflg(i) = .false.
 
  773 c  specify upper limit of mass flux at cloud base
 
  781           dp = 1000. * del(i,k)
 
  782           xmbmax(i) = dp / (g * dt2)
 
  789 c  compute cloud moisture property and precipitation
 
  795           qcko(i,kb(i)) = qo(i,kb(i))
 
  796           qrcko(i,kb(i)) = qo(i,kb(i))
 
  804             if(k.gt.kb(i).and.k.lt.ktcon(i)) 
then 
  805               dz    = zi(i,k) - zi(i,k-1)
 
  806               gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
  808      &             + gamma * dbyo(i,k) / (hvap * (1. + gamma))
 
  810               tem  = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
 
  811               tem1 = 0.5 * xlamud(i) * dz
 
  812               factor = 1. + tem - tem1
 
  813               qcko(i,k) = ((1.-tem1)*qcko(i,k-1)+tem*0.5*
 
  814      &                     (qo(i,k)+qo(i,k-1)))/factor
 
  815               qrcko(i,k) = qcko(i,k)
 
  817               dq = eta(i,k) * (qcko(i,k) - qrch)
 
  821 c  check 
if there is excess moisture to release latent heat
 
  823               if(k.ge.kbcon(i).and.dq.gt.0.) 
then 
  824                 etah = .5 * (eta(i,k) + eta(i,k-1))
 
  825                 if(ncloud.gt.0..and.k.gt.jmin(i)) 
then 
  826                   dp = 1000. * del(i,k)
 
  827                   qlk = dq / (eta(i,k) + etah * (c0 + c1) * dz)
 
  828                   dellal(i,k) = etah * c1 * dz * qlk * g / dp
 
  830                   qlk = dq / (eta(i,k) + etah * c0 * dz)
 
  832                 aa1(i) = aa1(i) - dz * g * qlk
 
  833                 qcko(i,k) = qlk + qrch
 
  834                 pwo(i,k) = etah * c0 * dz * qlk
 
  835                 pwavo(i) = pwavo(i) + pwo(i,k)
 
  837                 cnvwt(i,k) = etah * qlk * g / dp
 
  851 c  calculate cloud work function
 
  861             if(k.ge.kbcon(i).and.k.lt.ktcon(i)) 
then 
  862               dz1 = zo(i,k+1) - zo(i,k)
 
  863               gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
  864               rfact =  1. + delta * cp * gamma
 
  867      &                 dz1 * (g / (cp * to(i,k)))
 
  868      &                 * dbyo(i,k) / (1. + gamma)
 
  873      &                 max(val,(qeso(i,k) - qo(i,k)))
 
  880         if(cnvflg(i).and.aa1(i).le.0.) cnvflg(i) = .false.
 
  885         totflg = totflg .and. (.not. cnvflg(i))
 
  890 c  estimate the onvective overshooting as the level
 
  891 c    
where the [aafac * cloud work function] becomes zero,
 
  892 c    which is the final cloud top
 
  897           aa2(i) = aafac * aa1(i)
 
  903         ktcon1(i) = kmax(i) - 1
 
  908             if(k.ge.ktcon(i).and.k.lt.kmax(i)) 
then 
  909               dz1 = zo(i,k+1) - zo(i,k)
 
  910               gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
  911               rfact =  1. + delta * cp * gamma
 
  914      &                 dz1 * (g / (cp * to(i,k)))
 
  915      &                 * dbyo(i,k) / (1. + gamma)
 
  917               if(aa2(i).lt.0.) 
then 
  926 c  compute cloud moisture property, detraining cloud water
 
  927 c    and precipitation in overshooting layers
 
  933             if(k.ge.ktcon(i).and.k.lt.ktcon1(i)) 
then 
  934               dz    = zi(i,k) - zi(i,k-1)
 
  935               gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
  937      &             + gamma * dbyo(i,k) / (hvap * (1. + gamma))
 
  939               tem  = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
 
  940               tem1 = 0.5 * xlamud(i) * dz
 
  941               factor = 1. + tem - tem1
 
  942               qcko(i,k) = ((1.-tem1)*qcko(i,k-1)+tem*0.5*
 
  943      &                     (qo(i,k)+qo(i,k-1)))/factor
 
  944               qrcko(i,k) = qcko(i,k)
 
  946               dq = eta(i,k) * (qcko(i,k) - qrch)
 
  948 c  check 
if there is excess moisture to release latent heat
 
  951                 etah = .5 * (eta(i,k) + eta(i,k-1))
 
  952                 if(ncloud.gt.0.) 
then 
  953                   dp = 1000. * del(i,k)
 
  954                   qlk = dq / (eta(i,k) + etah * (c0 + c1) * dz)
 
  955                   dellal(i,k) = etah * c1 * dz * qlk * g / dp
 
  957                   qlk = dq / (eta(i,k) + etah * c0 * dz)
 
  959                 qcko(i,k) = qlk + qrch
 
  960                 pwo(i,k) = etah * c0 * dz * qlk
 
  961                 pwavo(i) = pwavo(i) + pwo(i,k)
 
  963                 cnvwt(i,k) = etah * qlk * g / dp
 
  970 c exchange ktcon with ktcon1
 
  981 c  this section is ready for cloud water
 
  986 c  compute liquid and vapor separation at cloud top
 
  991           gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
  993      &         + gamma * dbyo(i,k) / (hvap * (1. + gamma))
 
  994           dq = qcko(i,k) - qrch
 
  996 c  check 
if there is excess moisture to release latent heat
 
 1006 ccccc 
if(lat.eq.latd.and.lon.eq.lond.and.cnvflg(i)) 
then 
 1007 ccccc   print *, 
' aa1(i) before dwndrft =', aa1(i)
 
 1010 c------- downdraft calculations
 
 1012 c--- compute precipitation efficiency in terms of windshear
 
 1028             if(k.gt.kb(i).and.k.le.ktcon(i)) 
then 
 1029               shear= sqrt((uo(i,k)-uo(i,k-1)) ** 2
 
 1030      &                  + (vo(i,k)-vo(i,k-1)) ** 2)
 
 1031               vshear(i) = vshear(i) + shear
 
 1038           vshear(i) = 1.e3 * vshear(i) / (zi(i,ktcon(i))-zi(i,kb(i)))
 
 1039           e1=1.591-.639*vshear(i)
 
 1040      &       +.0953*(vshear(i)**2)-.00496*(vshear(i)**3)
 
 1043           edt(i) = min(edt(i),val)
 
 1045           edt(i) = max(edt(i),val)
 
 1051 c  determine detrainment rate between 1 and kbcon
 
 1065         if(cnvflg(i).and.k.ge.1.and.k.lt.kbcon(i)) 
then 
 1066           dz = zi(i,k+1) - zi(i,k)
 
 1067           sumx(i) = sumx(i) + dz
 
 1073         if(islimsk(i) == 1) beta = betal
 
 1075           dz  = (sumx(i)+zi(i,1))/float(kbcon(i))
 
 1076           tem = 1./float(kbcon(i))
 
 1077           xlamd(i) = (1.-beta**tem)/dz
 
 1081 c  determine downdraft mass flux
 
 1086           if (cnvflg(i) .and. k .le. kmax(i)-1) 
then 
 1087            if(k.lt.jmin(i).and.k.ge.kbcon(i)) 
then 
 1088               dz        = zi(i,k+1) - zi(i,k)
 
 1089               ptem      = xlamdd - xlamde
 
 1090               etad(i,k) = etad(i,k+1) * (1. - ptem * dz)
 
 1091            else if(k.lt.kbcon(i)) 
then 
 1092               dz        = zi(i,k+1) - zi(i,k)
 
 1093               ptem      = xlamd(i) + xlamdd - xlamde
 
 1094               etad(i,k) = etad(i,k+1) * (1. - ptem * dz)
 
 1100 c--- downdraft moisture properties
 
 1106           hcdo(i,jmn) = heo(i,jmn)
 
 1107           qcdo(i,jmn) = qo(i,jmn)
 
 1108           qrcdo(i,jmn)= qo(i,jmn)
 
 1109           ucdo(i,jmn) = uo(i,jmn)
 
 1110           vcdo(i,jmn) = vo(i,jmn)
 
 1118           if (cnvflg(i) .and. k.lt.jmin(i)) 
then 
 1119               dz = zi(i,k+1) - zi(i,k)
 
 1120               if(k.ge.kbcon(i)) 
then 
 1122                  tem1 = 0.5 * xlamdd * dz
 
 1125                  tem1 = 0.5 * (xlamd(i)+xlamdd) * dz
 
 1127               factor = 1. + tem - tem1
 
 1128               ptem = 0.5 * tem - pgcon
 
 1129               ptem1= 0.5 * tem + pgcon
 
 1130               hcdo(i,k) = ((1.-tem1)*hcdo(i,k+1)+tem*0.5*
 
 1131      &                     (heo(i,k)+heo(i,k+1)))/factor
 
 1132               ucdo(i,k) = ((1.-tem1)*ucdo(i,k+1)+ptem*uo(i,k+1)
 
 1133      &                     +ptem1*uo(i,k))/factor
 
 1134               vcdo(i,k) = ((1.-tem1)*vcdo(i,k+1)+ptem*vo(i,k+1)
 
 1135      &                     +ptem1*vo(i,k))/factor
 
 1136               dbyo(i,k) = hcdo(i,k) - heso(i,k)
 
 1144           if (cnvflg(i).and.k.lt.jmin(i)) 
then 
 1145               gamma      = el2orc * qeso(i,k) / (to(i,k)**2)
 
 1146               qrcdo(i,k) = qeso(i,k)+
 
 1147      &                (1./hvap)*(gamma/(1.+gamma))*dbyo(i,k)
 
 1150               dz = zi(i,k+1) - zi(i,k)
 
 1151               if(k.ge.kbcon(i)) 
then 
 1153                  tem1 = 0.5 * xlamdd * dz
 
 1156                  tem1 = 0.5 * (xlamd(i)+xlamdd) * dz
 
 1158               factor = 1. + tem - tem1
 
 1159               qcdo(i,k) = ((1.-tem1)*qrcdo(i,k+1)+tem*0.5*
 
 1160      &                     (qo(i,k)+qo(i,k+1)))/factor
 
 1167               pwdo(i,k)  = etad(i,k) * (qcdo(i,k) - qrcdo(i,k))
 
 1168               pwevo(i)   = pwevo(i) + pwdo(i,k)
 
 1173 c--- final downdraft strength dependent on precip
 
 1174 c--- efficiency(edt), normalized condensate(pwav), and
 
 1175 c--- evaporate(pwev)
 
 1180         if(islimsk(i) == 0) edtmax = edtmaxs
 
 1182           if(pwevo(i).lt.0.) 
then 
 1183             edto(i) = -edto(i) * pwavo(i) / pwevo(i)
 
 1184             edto(i) = min(edto(i),edtmax)
 
 1191 c--- downdraft cloudwork functions
 
 1196           if (cnvflg(i) .and. k .lt. jmin(i)) 
then 
 1197               gamma = el2orc * qeso(i,k) / to(i,k)**2
 
 1202               dz=-1.*(zo(i,k+1)-zo(i,k))
 
 1203               aa1(i)=aa1(i)+edto(i)*dz*(g/(cp*dt))*((dhh-dh)/(1.+dg))
 
 1204      &               *(1.+delta*cp*dg*dt/hvap)
 
 1206               aa1(i)=aa1(i)+edto(i)*
 
 1207      &        dz*g*delta*max(val,(qeso(i,k)-qo(i,k)))
 
 1213         if(cnvflg(i).and.aa1(i).le.0.) 
then 
 1220         totflg = totflg .and. (.not. cnvflg(i))
 
 1225 c--- what would the change be, that a cloud with unit mass
 
 1226 c--- will 
do to the environment?
 
 1231           if(cnvflg(i) .and. k .le. kmax(i)) 
then 
 1241           dp = 1000. * del(i,1)
 
 1242           dellah(i,1) = edto(i) * etad(i,1) * (hcdo(i,1)
 
 1243      &                   - heo(i,1)) * g / dp
 
 1244           dellaq(i,1) = edto(i) * etad(i,1) * (qrcdo(i,1)
 
 1245      &                   - qo(i,1)) * g / dp
 
 1246           dellau(i,1) = edto(i) * etad(i,1) * (ucdo(i,1)
 
 1247      &                   - uo(i,1)) * g / dp
 
 1248           dellav(i,1) = edto(i) * etad(i,1) * (vcdo(i,1)
 
 1249      &                   - vo(i,1)) * g / dp
 
 1253 c--- changed due to subsidence and entrainment
 
 1257           if (cnvflg(i).and.k.lt.ktcon(i)) 
then 
 1259               if(k.le.kb(i)) aup = 0.
 
 1261               if(k.gt.jmin(i)) adw = 0.
 
 1262               dp = 1000. * del(i,k)
 
 1263               dz = zi(i,k) - zi(i,k-1)
 
 1266               dv2h = .5 * (heo(i,k) + heo(i,k-1))
 
 1269               dv2q = .5 * (qo(i,k) + qo(i,k-1))
 
 1272               dv2u = .5 * (uo(i,k) + uo(i,k-1))
 
 1275               dv2v = .5 * (vo(i,k) + vo(i,k-1))
 
 1278               tem  = 0.5 * (xlamue(i,k)+xlamue(i,k-1))
 
 1281               if(k.le.kbcon(i)) 
then 
 1283                 ptem1 = xlamd(i)+xlamdd
 
 1289               dellah(i,k) = dellah(i,k) +
 
 1290      &     ((aup*eta(i,k)-adw*edto(i)*etad(i,k))*dv1h
 
 1291      &    - (aup*eta(i,k-1)-adw*edto(i)*etad(i,k-1))*dv3h
 
 1292      &    - (aup*tem*eta(i,k-1)+adw*edto(i)*ptem*etad(i,k))*dv2h*dz
 
 1293      &    +  aup*tem1*eta(i,k-1)*.5*(hcko(i,k)+hcko(i,k-1))*dz
 
 1294      &    +  adw*edto(i)*ptem1*etad(i,k)*.5*(hcdo(i,k)+hcdo(i,k-1))*dz
 
 1297               dellaq(i,k) = dellaq(i,k) +
 
 1298      &     ((aup*eta(i,k)-adw*edto(i)*etad(i,k))*dv1q
 
 1299      &    - (aup*eta(i,k-1)-adw*edto(i)*etad(i,k-1))*dv3q
 
 1300      &    - (aup*tem*eta(i,k-1)+adw*edto(i)*ptem*etad(i,k))*dv2q*dz
 
 1301      &    +  aup*tem1*eta(i,k-1)*.5*(qrcko(i,k)+qcko(i,k-1))*dz
 
 1302      &    +  adw*edto(i)*ptem1*etad(i,k)*.5*(qrcdo(i,k)+qcdo(i,k-1))*dz
 
 1305               dellau(i,k) = dellau(i,k) +
 
 1306      &     ((aup*eta(i,k)-adw*edto(i)*etad(i,k))*dv1u
 
 1307      &    - (aup*eta(i,k-1)-adw*edto(i)*etad(i,k-1))*dv3u
 
 1308      &    - (aup*tem*eta(i,k-1)+adw*edto(i)*ptem*etad(i,k))*dv2u*dz
 
 1309      &    +  aup*tem1*eta(i,k-1)*.5*(ucko(i,k)+ucko(i,k-1))*dz
 
 1310      &    +  adw*edto(i)*ptem1*etad(i,k)*.5*(ucdo(i,k)+ucdo(i,k-1))*dz
 
 1311      &    -  pgcon*(aup*eta(i,k-1)-adw*edto(i)*etad(i,k))*(dv1u-dv3u)
 
 1314               dellav(i,k) = dellav(i,k) +
 
 1315      &     ((aup*eta(i,k)-adw*edto(i)*etad(i,k))*dv1v
 
 1316      &    - (aup*eta(i,k-1)-adw*edto(i)*etad(i,k-1))*dv3v
 
 1317      &    - (aup*tem*eta(i,k-1)+adw*edto(i)*ptem*etad(i,k))*dv2v*dz
 
 1318      &    +  aup*tem1*eta(i,k-1)*.5*(vcko(i,k)+vcko(i,k-1))*dz
 
 1319      &    +  adw*edto(i)*ptem1*etad(i,k)*.5*(vcdo(i,k)+vcdo(i,k-1))*dz
 
 1320      &    -  pgcon*(aup*eta(i,k-1)-adw*edto(i)*etad(i,k))*(dv1v-dv3v)
 
 1332           dp = 1000. * del(i,indx)
 
 1333           dv1h = heo(i,indx-1)
 
 1334           dellah(i,indx) = eta(i,indx-1) *
 
 1335      &                     (hcko(i,indx-1) - dv1h) * g / dp
 
 1337           dellaq(i,indx) = eta(i,indx-1) *
 
 1338      &                     (qcko(i,indx-1) - dv1q) * g / dp
 
 1340           dellau(i,indx) = eta(i,indx-1) *
 
 1341      &                     (ucko(i,indx-1) - dv1u) * g / dp
 
 1343           dellav(i,indx) = eta(i,indx-1) *
 
 1344      &                     (vcko(i,indx-1) - dv1v) * g / dp
 
 1348           dellal(i,indx) = eta(i,indx-1) *
 
 1349      &                     qlko_ktcon(i) * g / dp
 
 1353 c------- final changed variable per unit mass flux
 
 1358           if (cnvflg(i).and.k .le. kmax(i)) 
then 
 1359             if(k.gt.ktcon(i)) 
then 
 1363             if(k.le.ktcon(i)) 
then 
 1364               qo(i,k) = dellaq(i,k) * mbdt + q1(i,k)
 
 1365               dellat = (dellah(i,k) - hvap * dellaq(i,k)) / cp
 
 1366               to(i,k) = dellat * mbdt + t1(i,k)
 
 1368               qo(i,k) = max(qo(i,k), val  )
 
 1375 c--- the above changed environment is now used to calulate the
 
 1376 c--- effect the arbitrary cloud(with unit mass flux)
 
 1377 c--- would have on the stability,
 
 1378 c--- which 
then is used to calculate the 
real mass flux,
 
 1379 c--- necessary to keep this change in balance with the large-scale
 
 1380 c--- destabilization.
 
 1382 c--- environmental conditions again, first heights
 
 1389           if(cnvflg(i) .and. k .le. kmax(i)) 
then 
 1390             qeso(i,k) = 0.01 * fpvs(to(i,k))      
 
 1391             qeso(i,k) = eps * qeso(i,k) / (pfld(i,k)+epsm1*qeso(i,k))
 
 1393             qeso(i,k) = max(qeso(i,k), val )
 
 1399 c--- moist static energy
 
 1404           if(cnvflg(i) .and. k .le. kmax(i)-1) 
then 
 1405             dz = .5 * (zo(i,k+1) - zo(i,k))
 
 1406             dp = .5 * (pfld(i,k+1) - pfld(i,k))
 
 1407             es = 0.01 * fpvs(to(i,k+1))      
 
 1408             pprime = pfld(i,k+1) + epsm1 * es
 
 1409             qs = eps * es / pprime
 
 1410             dqsdp = - qs / pprime
 
 1411             desdt = es * (fact1 / to(i,k+1) + fact2 / (to(i,k+1)**2))
 
 1412             dqsdt = qs * pfld(i,k+1) * desdt / (es * pprime)
 
 1413             gamma = el2orc * qeso(i,k+1) / (to(i,k+1)**2)
 
 1414             dt = (g * dz + hvap * dqsdp * dp) / (cp * (1. + gamma))
 
 1415             dq = dqsdt * dt + dqsdp * dp
 
 1416             to(i,k) = to(i,k+1) + dt
 
 1417             qo(i,k) = qo(i,k+1) + dq
 
 1418             po(i,k) = .5 * (pfld(i,k) + pfld(i,k+1))
 
 1424           if(cnvflg(i) .and. k .le. kmax(i)-1) 
then 
 1425             qeso(i,k) = 0.01 * fpvs(to(i,k))      
 
 1426             qeso(i,k) = eps * qeso(i,k) / (po(i,k) + epsm1 * qeso(i,k))
 
 1428             qeso(i,k) = max(qeso(i,k), val1)
 
 1430             qo(i,k)   = max(qo(i,k), val2 )
 
 1432             heo(i,k)   = .5 * g * (zo(i,k) + zo(i,k+1)) +
 
 1433      &                    cp * to(i,k) + hvap * qo(i,k)
 
 1434             heso(i,k) = .5 * g * (zo(i,k) + zo(i,k+1)) +
 
 1435      &                  cp * to(i,k) + hvap * qeso(i,k)
 
 1442           heo(i,k) = g * zo(i,k) + cp * to(i,k) + hvap * qo(i,k)
 
 1443           heso(i,k) = g * zo(i,k) + cp * to(i,k) + hvap * qeso(i,k)
 
 1444 c         heo(i,k) = min(heo(i,k),heso(i,k))
 
 1448 c**************************** static control
 
 1450 c------- moisture and cloud work functions
 
 1463           hcko(i,indx) = heo(i,indx)
 
 1464           qcko(i,indx) = qo(i,indx)
 
 1470             if(k.gt.kb(i).and.k.le.ktcon(i)) 
then 
 1471               dz = zi(i,k) - zi(i,k-1)
 
 1472               tem  = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
 
 1473               tem1 = 0.5 * xlamud(i) * dz
 
 1474               factor = 1. + tem - tem1
 
 1475               hcko(i,k) = ((1.-tem1)*hcko(i,k-1)+tem*0.5*
 
 1476      &                     (heo(i,k)+heo(i,k-1)))/factor
 
 1484             if(k.gt.kb(i).and.k.lt.ktcon(i)) 
then 
 1485               dz = zi(i,k) - zi(i,k-1)
 
 1486               gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
 1487               xdby = hcko(i,k) - heso(i,k)
 
 1489      &              + gamma * xdby / (hvap * (1. + gamma))
 
 1491               tem  = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
 
 1492               tem1 = 0.5 * xlamud(i) * dz
 
 1493               factor = 1. + tem - tem1
 
 1494               qcko(i,k) = ((1.-tem1)*qcko(i,k-1)+tem*0.5*
 
 1495      &                     (qo(i,k)+qo(i,k-1)))/factor
 
 1497               dq = eta(i,k) * (qcko(i,k) - xqrch)
 
 1499               if(k.ge.kbcon(i).and.dq.gt.0.) 
then 
 1500                 etah = .5 * (eta(i,k) + eta(i,k-1))
 
 1501                 if(ncloud.gt.0..and.k.gt.jmin(i)) 
then 
 1502                   qlk = dq / (eta(i,k) + etah * (c0 + c1) * dz)
 
 1504                   qlk = dq / (eta(i,k) + etah * c0 * dz)
 
 1506                 if(k.lt.ktcon1(i)) 
then 
 1507                   xaa0(i) = xaa0(i) - dz * g * qlk
 
 1509                 qcko(i,k) = qlk + xqrch
 
 1510                 xpw = etah * c0 * dz * qlk
 
 1511                 xpwav(i) = xpwav(i) + xpw
 
 1514             if(k.ge.kbcon(i).and.k.lt.ktcon1(i)) 
then 
 1515               dz1 = zo(i,k+1) - zo(i,k)
 
 1516               gamma = el2orc * qeso(i,k) / (to(i,k)**2)
 
 1517               rfact =  1. + delta * cp * gamma
 
 1520      &                + dz1 * (g / (cp * to(i,k)))
 
 1521      &                * xdby / (1. + gamma)
 
 1526      &                 max(val,(qeso(i,k) - qo(i,k)))
 
 1532 c------- downdraft calculations
 
 1534 c--- downdraft moisture properties
 
 1540           hcdo(i,jmn) = heo(i,jmn)
 
 1541           qcdo(i,jmn) = qo(i,jmn)
 
 1542           qrcd(i,jmn) = qo(i,jmn)
 
 1549           if (cnvflg(i) .and. k.lt.jmin(i)) 
then 
 1550               dz = zi(i,k+1) - zi(i,k)
 
 1551               if(k.ge.kbcon(i)) 
then 
 1553                  tem1 = 0.5 * xlamdd * dz
 
 1556                  tem1 = 0.5 * (xlamd(i)+xlamdd) * dz
 
 1558               factor = 1. + tem - tem1
 
 1559               hcdo(i,k) = ((1.-tem1)*hcdo(i,k+1)+tem*0.5*
 
 1560      &                     (heo(i,k)+heo(i,k+1)))/factor
 
 1567           if (cnvflg(i) .and. k .lt. jmin(i)) 
then 
 1570               gamma    = el2orc * dq / dt**2
 
 1571               dh       = hcdo(i,k) - heso(i,k)
 
 1572               qrcd(i,k)=dq+(1./hvap)*(gamma/(1.+gamma))*dh
 
 1575               dz = zi(i,k+1) - zi(i,k)
 
 1576               if(k.ge.kbcon(i)) 
then 
 1578                  tem1 = 0.5 * xlamdd * dz
 
 1581                  tem1 = 0.5 * (xlamd(i)+xlamdd) * dz
 
 1583               factor = 1. + tem - tem1
 
 1584               qcdo(i,k) = ((1.-tem1)*qrcd(i,k+1)+tem*0.5*
 
 1585      &                     (qo(i,k)+qo(i,k+1)))/factor
 
 1592               xpwd     = etad(i,k) * (qcdo(i,k) - qrcd(i,k))
 
 1593               xpwev(i) = xpwev(i) + xpwd
 
 1600         if(islimsk(i) == 0) edtmax = edtmaxs
 
 1602           if(xpwev(i).ge.0.) 
then 
 1605             edtx(i) = -edtx(i) * xpwav(i) / xpwev(i)
 
 1606             edtx(i) = min(edtx(i),edtmax)
 
 1612 c--- downdraft cloudwork functions
 
 1617           if (cnvflg(i) .and. k.lt.jmin(i)) 
then 
 1618               gamma = el2orc * qeso(i,k) / to(i,k)**2
 
 1623               dz=-1.*(zo(i,k+1)-zo(i,k))
 
 1624               xaa0(i)=xaa0(i)+edtx(i)*dz*(g/(cp*dt))*((dhh-dh)/(1.+dg))
 
 1625      &                *(1.+delta*cp*dg*dt/hvap)
 
 1627               xaa0(i)=xaa0(i)+edtx(i)*
 
 1628      &        dz*g*delta*max(val,(qeso(i,k)-qo(i,k)))
 
 1633 c  calculate critical cloud work function
 
 1639           if(pfld(i,ktcon(i)).lt.pcrit(15))
then 
 1640             acrt(i)=acrit(15)*(975.-pfld(i,ktcon(i)))
 
 1642           else if(pfld(i,ktcon(i)).gt.pcrit(1))
then 
 1645             k =  int((850. - pfld(i,ktcon(i)))/50.) + 2
 
 1648             acrt(i)=acrit(k)+(acrit(k-1)-acrit(k))*
 
 1649      &           (pfld(i,ktcon(i))-pcrit(k))/(pcrit(k-1)-pcrit(k))
 
 1656           if(islimsk(i) == 1) 
then 
 1668 c  modify critical cloud workfunction by cloud base vertical velocity
 
 1670           if(pdot(i).le.w4) 
then 
 1671             acrtfct(i) = (pdot(i) - w4) / (w3 - w4)
 
 1672           elseif(pdot(i).ge.-w4) 
then 
 1673             acrtfct(i) = - (pdot(i) + w4) / (w4 - w3)
 
 1678           acrtfct(i) = max(acrtfct(i),val1)
 
 1680           acrtfct(i) = min(acrtfct(i),val2)
 
 1681           acrtfct(i) = 1. - acrtfct(i)
 
 1683 c  modify acrtfct(i) by colume mean rh 
if rhbar(i) is greater than 80 percent
 
 1685 c         
if(rhbar(i).ge..8) 
then 
 1686 c           acrtfct(i) = acrtfct(i) * (.9 - min(rhbar(i),.9)) * 10.
 
 1689 c  modify adjustment time scale by cloud base vertical velocity
 
 1692           dtconv(i) = dt2 + max((1800. - dt2),0.) *
 
 1693      &                (pdot(i) - w2) / (w1 - w2)
 
 1694 c         dtconv(i) = max(dtconv(i), dt2)
 
 1695 c         dtconv(i) = 1800. * (pdot(i) - w2) / (w1 - w2)
 
 1696           dtconv(i) = max(dtconv(i),dtmin)
 
 1697           dtconv(i) = min(dtconv(i),dtmax)
 
 1702 c--- large scale forcing
 
 1711           fld(i)=(aa1(i)-acrt(i)* acrtfct(i))/dtconv(i)
 
 1712           if(fld(i).le.0.) cnvflg(i) = .false.
 
 1720 c         xaa0(i) = max(xaa0(i),0.)
 
 1721           xk(i) = (xaa0(i) - aa1(i)) / mbdt
 
 1722           if(xk(i).ge.0.) cnvflg(i) = .false.
 
 1725 c--- kernel, cloud base mass flux
 
 1732           xmb(i) = -fld(i) / xk(i)
 
 1733           xmb(i) = min(xmb(i),xmbmax(i))
 
 1740         totflg = totflg .and. (.not. cnvflg(i))
 
 1745 c  restore to,qo,uo,vo to t1,q1,u1,v1 in 
case convection stops
 
 1750           if (cnvflg(i) .and. k .le. kmax(i)) 
then 
 1755             qeso(i,k) = 0.01 * fpvs(t1(i,k))      
 
 1756             qeso(i,k) = eps * qeso(i,k) / (pfld(i,k) + epsm1*qeso(i,k))
 
 1758             qeso(i,k) = max(qeso(i,k), val )
 
 1764 c--- feedback: simply the changes from the cloud with unit mass flux
 
 1765 c---           multiplied by  the mass flux necessary to keep the
 
 1766 c---           equilibrium with the larger-scale.
 
 1782           if (cnvflg(i) .and. k .le. kmax(i)) 
then 
 1783             if(k.le.ktcon(i)) 
then 
 1784               dellat = (dellah(i,k) - hvap * dellaq(i,k)) / cp
 
 1785               t1(i,k) = t1(i,k) + dellat * xmb(i) * dt2
 
 1786               q1(i,k) = q1(i,k) + dellaq(i,k) * xmb(i) * dt2
 
 1790               u1(i,k) = u1(i,k) + dellau(i,k) * xmb(i) * dt2
 
 1791               v1(i,k) = v1(i,k) + dellav(i,k) * xmb(i) * dt2
 
 1792               dp = 1000. * del(i,k)
 
 1793               delhbar(i) = delhbar(i) + dellah(i,k)*xmb(i)*dp/g
 
 1794               delqbar(i) = delqbar(i) + dellaq(i,k)*xmb(i)*dp/g
 
 1795               deltbar(i) = deltbar(i) + dellat*xmb(i)*dp/g
 
 1796               delubar(i) = delubar(i) + dellau(i,k)*xmb(i)*dp/g
 
 1797               delvbar(i) = delvbar(i) + dellav(i,k)*xmb(i)*dp/g
 
 1805           if (cnvflg(i) .and. k .le. kmax(i)) 
then 
 1806             if(k.le.ktcon(i)) 
then 
 1807               qeso(i,k) = 0.01 * fpvs(t1(i,k))      
 
 1808               qeso(i,k) = eps * qeso(i,k)/(pfld(i,k) + epsm1*qeso(i,k))
 
 1810               qeso(i,k) = max(qeso(i,k), val )
 
 1825           if (cnvflg(i) .and. k .le. kmax(i)) 
then 
 1826             if(k.lt.ktcon(i)) 
then 
 1828               if(k.le.kb(i)) aup = 0.
 
 1830               if(k.ge.jmin(i)) adw = 0.
 
 1831               rain =  aup * pwo(i,k) + adw * edto(i) * pwdo(i,k)
 
 1832               rntot(i) = rntot(i) + rain * xmb(i) * .001 * dt2
 
 1842           if (k .le. kmax(i)) 
then 
 1846             if(cnvflg(i).and.k.lt.ktcon(i)) 
then 
 1848               if(k.le.kb(i)) aup = 0.
 
 1850               if(k.ge.jmin(i)) adw = 0.
 
 1851               rain =  aup * pwo(i,k) + adw * edto(i) * pwdo(i,k)
 
 1852               rn(i) = rn(i) + rain * xmb(i) * .001 * dt2
 
 1854             if(flg(i).and.k.lt.ktcon(i)) 
then 
 1855               evef = edt(i) * evfact
 
 1856               if(islimsk(i) == 1) evef=edt(i) * evfactl
 
 1858 c             
if(islimsk(i) == 1) evef = 0.
 
 1859               qcond(i) = evef * (q1(i,k) - qeso(i,k))
 
 1860      &                 / (1. + el2orc * qeso(i,k) / t1(i,k)**2)
 
 1861               dp = 1000. * del(i,k)
 
 1862               if(rn(i).gt.0..and.qcond(i).lt.0.) 
then 
 1863                 qevap(i) = -qcond(i) * (1.-exp(-.32*sqrt(dt2*rn(i))))
 
 1864                 qevap(i) = min(qevap(i), rn(i)*1000.*g/dp)
 
 1865                 delq2(i) = delqev(i) + .001 * qevap(i) * dp / g
 
 1867               if(rn(i).gt.0..and.qcond(i).lt.0..and.
 
 1868      &           delq2(i).gt.rntot(i)) 
then 
 1869                 qevap(i) = 1000.* g * (rntot(i) - delqev(i)) / dp
 
 1872               if(rn(i).gt.0..and.qevap(i).gt.0.) 
then 
 1873                 q1(i,k) = q1(i,k) + qevap(i)
 
 1874                 t1(i,k) = t1(i,k) - elocp * qevap(i)
 
 1875                 rn(i) = rn(i) - .001 * qevap(i) * dp / g
 
 1876                 deltv(i) = - elocp*qevap(i)/dt2
 
 1877                 delq(i) =  + qevap(i)/dt2
 
 1878                 delqev(i) = delqev(i) + .001*dp*qevap(i)/g
 
 1880               dellaq(i,k) = dellaq(i,k) + delq(i) / xmb(i)
 
 1881               delqbar(i) = delqbar(i) + delq(i)*dp/g
 
 1882               deltbar(i) = deltbar(i) + deltv(i)*dp/g
 
 1899 c  precipitation rate converted to actual precip
 
 1900 c  in unit of m instead of kg
 
 1905 c  in the event of upper level rain evaporation and lower level downdraft
 
 1906 c    moistening, rn can become negative, in this 
case, we back out of the
 
 1907 c    heating and the moistening
 
 1910           if(rn(i).lt.0..and..not.flg(i)) rn(i) = 0.
 
 1911           if(rn(i).le.0.) 
then 
 1922 c  convective cloud water
 
 1927           if (cnvflg(i) .and. rn(i).gt.0.) 
then 
 1928             if (k.ge.kbcon(i).and.k.lt.ktcon(i)) 
then 
 1929               cnvw(i,k) = cnvwt(i,k) * xmb(i) * dt2
 
 1935 c  convective cloud cover
 
 1940           if (cnvflg(i) .and. rn(i).gt.0.) 
then 
 1941             if (k.ge.kbcon(i).and.k.lt.ktcon(i)) 
then 
 1942               cnvc(i,k) = 0.04 * log(1. + 675. * eta(i,k) * xmb(i))
 
 1943               cnvc(i,k) = min(cnvc(i,k), 0.6)
 
 1944               cnvc(i,k) = max(cnvc(i,k), 0.0)
 
 1954       if (ncloud.gt.0) 
then 
 1958           if (cnvflg(i) .and. rn(i).gt.0.) 
then 
 1959             if (k.gt.kb(i).and.k.le.ktcon(i)) 
then 
 1960               tem  = dellal(i,k) * xmb(i) * dt2
 
 1961               tem1 = max(0.0, min(1.0, (tcr-t1(i,k))*tcrf))
 
 1962               if (ql(i,k,2) .gt. -999.0) 
then 
 1963                 ql(i,k,1) = ql(i,k,1) + tem * tem1            
 
 1964                 ql(i,k,2) = ql(i,k,2) + tem *(1.0-tem1)       
 
 1966                 ql(i,k,1) = ql(i,k,1) + tem
 
 1978           if(cnvflg(i).and.rn(i).le.0.) 
then 
 1979             if (k .le. kmax(i)) 
then 
 1994           if(cnvflg(i).and.rn(i).gt.0.) 
then 
 1995             if(k.ge.kb(i) .and. k.lt.ktop(i)) 
then 
 1996               ud_mf(i,k) = eta(i,k) * xmb(i) * dt2
 
 2003         if(cnvflg(i).and.rn(i).gt.0.) 
then 
 2005            dt_mf(i,k) = ud_mf(i,k)
 
 2011           if(cnvflg(i).and.rn(i).gt.0.) 
then 
 2012             if(k.ge.1 .and. k.le.jmin(i)) 
then 
 2013               dd_mf(i,k) = edto(i) * etad(i,k) * xmb(i) * dt2
 
real(kind=kind_phys), parameter con_cliq
spec heat H2O liq ( ) 
 
subroutine sascnvn(im, ix, km, jcap, delt, delp, prslp, psp, phil, ql,                               q1, t1, u1, v1, cldwrk, rn, kbot, ktop, kcnv, islimsk,                                                                           dot, ncloud, ud_mf, dd_mf, dt_mf, cnvw, cnvc)
This subroutine contains the entirety of the SAS deep convection scheme. 
 
real(kind=kind_phys), parameter con_g
gravity ( ) 
 
real(kind=kind_phys), parameter con_rv
gas constant H2O ( ) 
 
real(kind=kind_phys), parameter con_t0c
temp at 0C (K) 
 
real(kind=kind_phys), parameter con_cvap
spec heat H2O gas ( ) 
 
real(kind=kind_phys), parameter con_cp
spec heat air at p ( ) 
 
real(kind=kind_phys), parameter con_hvap
lat heat H2O cond ( )