55 subroutine shalcnv(im,ix,km,jcap,delt,delp,prslp,psp,phil,ql, &
56 & q1,t1,u1,v1,rn,kbot,ktop,kcnv,islimsk, &
57 & dot,ncloud,hpbl,heat,evap,ud_mf,dt_mf,cnvw,cnvc)
61 use machine
, only : kind_phys
62 use funcphys
, only : fpvs
66 &, eps => con_eps, epsm1 => con_epsm1
69 integer im, ix, km, jcap, ncloud, &
70 & kbot(im), ktop(im), kcnv(im)
72 real(kind=kind_phys) delt
73 real(kind=kind_phys) psp(im), delp(ix,km), prslp(ix,km)
74 real(kind=kind_phys) ps(im), del(ix,km), prsl(ix,km), &
75 & ql(ix,km,2),q1(ix,km), t1(ix,km), &
76 & u1(ix,km), v1(ix,km), &
78 & dot(ix,km), phil(ix,km), hpbl(im), &
79 & heat(im), evap(im), cnvw(ix,km), &
81 & ud_mf(im,km),dt_mf(im,km)
83 integer i,j,indx, k, kk, km1
85 integer,
dimension(im),
intent(in) :: islimsk
87 real(kind=kind_phys) c0, dellat, delta,
90 & dq, dqsdp, dqsdt, dt,
91 & dt2, dtmax, dtmin, dv1h,
92 & dv1q, dv2h, dv2q, dv1u,
93 & dv1v, dv2u, dv2v, dv3q,
94 & dv3h, dv3u, dv3v, clam,
96 & el2orc, elocp, aafac,
98 & evef, evfact, evfactl, fact1,
99 & fact2, factor, fjcap,
100 & g, gamma, pprime, betaw,
102 & rfact, shear, tem1,
104 & val2, w1, w1l, w1s,
107 & w4s, tem, ptem, ptem1,
110 integer kb(im), kbcon(im), kbcon1(im),
111 & ktcon(im), ktcon1(im),
114 real(kind=kind_phys) aa1(im),
115 & delhbar(im), delq(im), delq2(im),
116 & delqbar(im), delqev(im), deltbar(im),
117 & deltv(im), edt(im),
118 & wstar(im), sflx(im),
119 & pdot(im), po(im,km),
120 & qcond(im), qevap(im), hmax(im),
121 & rntot(im), vshear(im),
122 & xlamud(im), xmb(im), xmbmax(im),
123 & delubar(im), delvbar(im)
125 real(kind=kind_phys) cincr, cincrmax, cincrmin
127 c physical parameters
129 parameter(elocp=hvap/cp,
130 & el2orc=hvap*hvap/(rv*cp))
131 parameter(c0=.002,c1=5.e-4,delta=fv)
132 parameter(fact1=(cvap-cliq)/rv,fact2=hvap/rv-fact1*t0c)
133 parameter(cincrmax=180.,cincrmin=120.,dthk=25.)
134 parameter(h1=0.33333333)
135 c local variables and arrays
136 real(kind=kind_phys) pfld(im,km), to(im,km), qo(im,km),
137 & uo(im,km), vo(im,km), qeso(im,km)
140 real(kind=kind_phys) qlko_ktcon(im), dellal(im,km),
141 & dbyo(im,km), zo(im,km), xlamue(im,km),
142 & heo(im,km), heso(im,km),
143 & dellah(im,km), dellaq(im,km),
144 & dellau(im,km), dellav(im,km), hcko(im,km),
145 & ucko(im,km), vcko(im,km), qcko(im,km),
146 & qrcko(im,km), eta(im,km),
147 & zi(im,km), pwo(im,km),
148 & tx1(im), cnvwt(im,km)
150 logical totflg, cnvflg(im), flg(im)
152 real(kind=kind_phys) tf, tcr, tcrf
153 parameter(tf=233.16, tcr=263.16, tcrf=1.0/(tcr-tf))
155 c-----------------------------------------------------------------------
168 c compute surface buoyancy flux
176 sflx(i) = heat(i)+fv*t1(i,1)*evap(i)
184 if(kcnv(i).eq.1) cnvflg(i) = .false.
185 if(sflx(i).le.0.) cnvflg(i) = .false.
212 totflg = totflg .and. (.not. cnvflg(i))
220 dtmin = max(dt2, val )
222 dtmax = max(dt2, val )
223 c model tunable parameters are all here
233 fjcap = (float(jcap) / 126.) ** 2
235 fjcap = max(fjcap,val)
245 c define top layer for search of the downdraft originating layer
246 c and the maximum thetae for updraft
257 if (prsl(i,k)*tx1(i) .gt. 0.70) kbm(i) = k + 1
258 if (prsl(i,k)*tx1(i) .gt. 0.60) kmax(i) = k + 1
262 kbm(i) = min(kbm(i),kmax(i))
265 c hydrostatic height assume zero terr and compute
266 c updraft entrainment rate as an inverse
function of height
271 zo(i,k) = phil(i,k) / g
277 zi(i,k) = 0.5*(zo(i,k)+zo(i,k+1))
278 xlamue(i,k) = clam / zi(i,k)
282 xlamue(i,km) = xlamue(i,km1)
294 if (flg(i).and.zo(i,k).le.hpbl(i))
then
302 kpbl(i)= min(kpbl(i),kbm(i))
306 c convert surface pressure to mb from cb
311 if (cnvflg(i) .and. k .le. kmax(i))
then
312 pfld(i,k) = prsl(i,k) * 10.0
334 c p is pressure of the layer(mb)
335 c t is temperature at t-dt(k)..tn
336 c q is mixing ratio at t-dt(kg/kg)..qn
337 c to is temperature at t+dt(k)... this is after advection and turbulan
338 c qo is mixing ratio at t+dt(kg/kg)..q1
343 if (cnvflg(i) .and. k .le. kmax(i))
then
344 qeso(i,k) = 0.01 * fpvs(to(i,k))
345 qeso(i,k) = eps * qeso(i,k) / (pfld(i,k) + epsm1*qeso(i,k))
347 qeso(i,k) = max(qeso(i,k), val1)
349 qo(i,k) = max(qo(i,k), val2 )
356 c compute moist static energy
361 if (cnvflg(i) .and. k .le. kmax(i))
then
363 tem = phil(i,k) + cp * to(i,k)
364 heo(i,k) = tem + hvap * qo(i,k)
365 heso(i,k) = tem + hvap * qeso(i,k)
366 c heo(i,k) = min(heo(i,k),heso(i,k))
371 c determine level with largest moist static energy within pbl
372 c this is the level
where updraft starts
384 if (cnvflg(i).and.k.le.kpbl(i))
then
385 if(heo(i,k).gt.hmax(i))
then
396 if (cnvflg(i) .and. k .le. kmax(i)-1)
then
397 dz = .5 * (zo(i,k+1) - zo(i,k))
398 dp = .5 * (pfld(i,k+1) - pfld(i,k))
399 es = 0.01 * fpvs(to(i,k+1))
400 pprime = pfld(i,k+1) + epsm1 * es
401 qs = eps * es / pprime
402 dqsdp = - qs / pprime
403 desdt = es * (fact1 / to(i,k+1) + fact2 / (to(i,k+1)**2))
404 dqsdt = qs * pfld(i,k+1) * desdt / (es * pprime)
405 gamma = el2orc * qeso(i,k+1) / (to(i,k+1)**2)
406 dt = (g * dz + hvap * dqsdp * dp) / (cp * (1. + gamma))
407 dq = dqsdt * dt + dqsdp * dp
408 to(i,k) = to(i,k+1) + dt
409 qo(i,k) = qo(i,k+1) + dq
410 po(i,k) = .5 * (pfld(i,k) + pfld(i,k+1))
418 if (cnvflg(i) .and. k .le. kmax(i)-1)
then
419 qeso(i,k) = 0.01 * fpvs(to(i,k))
420 qeso(i,k) = eps * qeso(i,k) / (po(i,k) + epsm1*qeso(i,k))
422 qeso(i,k) = max(qeso(i,k), val1)
424 qo(i,k) = max(qo(i,k), val2 )
426 heo(i,k) = .5 * g * (zo(i,k) + zo(i,k+1)) +
427 & cp * to(i,k) + hvap * qo(i,k)
428 heso(i,k) = .5 * g * (zo(i,k) + zo(i,k+1)) +
429 & cp * to(i,k) + hvap * qeso(i,k)
430 uo(i,k) = .5 * (uo(i,k) + uo(i,k+1))
431 vo(i,k) = .5 * (vo(i,k) + vo(i,k+1))
436 c look for the level of free convection as cloud base
440 if(flg(i)) kbcon(i) = kmax(i)
444 if (flg(i).and.k.lt.kbm(i))
then
445 if(k.gt.kb(i).and.heo(i,kb(i)).gt.heso(i,k))
then
456 if(kbcon(i).eq.kmax(i)) cnvflg(i) = .false.
462 totflg = totflg .and. (.not. cnvflg(i))
467 c determine critical convective inhibition
468 c as a
function of vertical velocity at cloud base.
474 pdot(i) = 0.01 * dot(i,kbcon(i))
479 if(islimsk(i) == 1)
then
490 if(pdot(i).le.w4)
then
491 ptem = (pdot(i) - w4) / (w3 - w4)
492 elseif(pdot(i).ge.-w4)
then
493 ptem = - (pdot(i) + w4) / (w4 - w3)
498 ptem = max(ptem,val1)
500 ptem = min(ptem,val2)
502 ptem1= .5*(cincrmax-cincrmin)
503 cincr = cincrmax - ptem * ptem1
504 tem1 = pfld(i,kb(i)) - pfld(i,kbcon(i))
505 if(tem1.gt.cincr)
then
513 totflg = totflg .and. (.not. cnvflg(i))
518 c assume the detrainment rate for the updrafts to be same as
519 c the entrainment rate at cloud base
524 xlamud(i) = xlamue(i,kbcon(i))
528 c determine updraft mass flux for the subcloud layers
538 if(k.lt.kbcon(i).and.k.ge.kb(i))
then
539 dz = zi(i,k+1) - zi(i,k)
540 ptem = 0.5*(xlamue(i,k)+xlamue(i,k+1))-xlamud(i)
541 eta(i,k) = eta(i,k+1) / (1. + ptem * dz)
547 c compute mass flux above cloud base
552 if(k.gt.kbcon(i).and.k.lt.kmax(i))
then
553 dz = zi(i,k) - zi(i,k-1)
554 ptem = 0.5*(xlamue(i,k)+xlamue(i,k-1))-xlamud(i)
555 eta(i,k) = eta(i,k-1) * (1 + ptem * dz)
561 c compute updraft cloud property
567 hcko(i,indx) = heo(i,indx)
568 ucko(i,indx) = uo(i,indx)
569 vcko(i,indx) = vo(i,indx)
577 if(k.gt.kb(i).and.k.lt.kmax(i))
then
578 dz = zi(i,k) - zi(i,k-1)
579 tem = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
580 tem1 = 0.5 * xlamud(i) * dz
581 factor = 1. + tem - tem1
582 ptem = 0.5 * tem + pgcon
583 ptem1= 0.5 * tem - pgcon
584 hcko(i,k) = ((1.-tem1)*hcko(i,k-1)+tem*0.5*
585 & (heo(i,k)+heo(i,k-1)))/factor
586 ucko(i,k) = ((1.-tem1)*ucko(i,k-1)+ptem*uo(i,k)
587 & +ptem1*uo(i,k-1))/factor
588 vcko(i,k) = ((1.-tem1)*vcko(i,k-1)+ptem*vo(i,k)
589 & +ptem1*vo(i,k-1))/factor
590 dbyo(i,k) = hcko(i,k) - heso(i,k)
596 c taking account into convection inhibition due to existence of
597 c dry layers below cloud base
606 if (flg(i).and.k.lt.kbm(i))
then
607 if(k.ge.kbcon(i).and.dbyo(i,k).gt.0.)
then
616 if(kbcon1(i).eq.kmax(i)) cnvflg(i) = .false.
621 tem = pfld(i,kbcon(i)) - pfld(i,kbcon1(i))
630 totflg = totflg .and. (.not. cnvflg(i))
635 c determine first guess cloud top as the level of zero buoyancy
636 c limited to the level of sigma=0.7
641 if(flg(i)) ktcon(i) = kbm(i)
645 if (flg(i).and.k .lt. kbm(i))
then
646 if(k.gt.kbcon1(i).and.dbyo(i,k).lt.0.)
then
654 c turn off shallow convection
if cloud top is less than pbl top
670 c specify upper limit of mass flux at cloud base
678 dp = 1000. * del(i,k)
679 xmbmax(i) = dp / (g * dt2)
686 c compute cloud moisture property and precipitation
692 qcko(i,kb(i)) = qo(i,kb(i))
693 qrcko(i,kb(i)) = qo(i,kb(i))
700 if(k.gt.kb(i).and.k.lt.ktcon(i))
then
701 dz = zi(i,k) - zi(i,k-1)
702 gamma = el2orc * qeso(i,k) / (to(i,k)**2)
704 & + gamma * dbyo(i,k) / (hvap * (1. + gamma))
706 tem = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
707 tem1 = 0.5 * xlamud(i) * dz
708 factor = 1. + tem - tem1
709 qcko(i,k) = ((1.-tem1)*qcko(i,k-1)+tem*0.5*
710 & (qo(i,k)+qo(i,k-1)))/factor
711 qrcko(i,k) = qcko(i,k)
713 dq = eta(i,k) * (qcko(i,k) - qrch)
717 c below lfc check
if there is excess moisture to release latent heat
719 if(k.ge.kbcon(i).and.dq.gt.0.)
then
720 etah = .5 * (eta(i,k) + eta(i,k-1))
721 if(ncloud.gt.0.)
then
722 dp = 1000. * del(i,k)
723 qlk = dq / (eta(i,k) + etah * (c0 + c1) * dz)
724 dellal(i,k) = etah * c1 * dz * qlk * g / dp
726 qlk = dq / (eta(i,k) + etah * c0 * dz)
728 aa1(i) = aa1(i) - dz * g * qlk
729 qcko(i,k)= qlk + qrch
730 pwo(i,k) = etah * c0 * dz * qlk
731 cnvwt(i,k) = etah * qlk * g / dp
738 c calculate cloud work function
748 if(k.ge.kbcon(i).and.k.lt.ktcon(i))
then
749 dz1 = zo(i,k+1) - zo(i,k)
750 gamma = el2orc * qeso(i,k) / (to(i,k)**2)
751 rfact = 1. + delta * cp * gamma
754 & dz1 * (g / (cp * to(i,k)))
755 & * dbyo(i,k) / (1. + gamma)
760 & max(val,(qeso(i,k) - qo(i,k)))
767 if(cnvflg(i).and.aa1(i).le.0.) cnvflg(i) = .false.
772 totflg = totflg .and. (.not. cnvflg(i))
777 c estimate the onvective overshooting as the level
778 c
where the [aafac * cloud work function] becomes zero,
779 c which is the final cloud top
780 c limited to the level of sigma=0.7
785 aa1(i) = aafac * aa1(i)
796 if(k.ge.ktcon(i).and.k.lt.kbm(i))
then
797 dz1 = zo(i,k+1) - zo(i,k)
798 gamma = el2orc * qeso(i,k) / (to(i,k)**2)
799 rfact = 1. + delta * cp * gamma
802 & dz1 * (g / (cp * to(i,k)))
803 & * dbyo(i,k) / (1. + gamma)
805 if(aa1(i).lt.0.)
then
814 c compute cloud moisture property, detraining cloud water
815 c and precipitation in overshooting layers
821 if(k.ge.ktcon(i).and.k.lt.ktcon1(i))
then
822 dz = zi(i,k) - zi(i,k-1)
823 gamma = el2orc * qeso(i,k) / (to(i,k)**2)
825 & + gamma * dbyo(i,k) / (hvap * (1. + gamma))
827 tem = 0.5 * (xlamue(i,k)+xlamue(i,k-1)) * dz
828 tem1 = 0.5 * xlamud(i) * dz
829 factor = 1. + tem - tem1
830 qcko(i,k) = ((1.-tem1)*qcko(i,k-1)+tem*0.5*
831 & (qo(i,k)+qo(i,k-1)))/factor
832 qrcko(i,k) = qcko(i,k)
834 dq = eta(i,k) * (qcko(i,k) - qrch)
836 c check
if there is excess moisture to release latent heat
839 etah = .5 * (eta(i,k) + eta(i,k-1))
840 if(ncloud.gt.0.)
then
841 dp = 1000. * del(i,k)
842 qlk = dq / (eta(i,k) + etah * (c0 + c1) * dz)
843 dellal(i,k) = etah * c1 * dz * qlk * g / dp
845 qlk = dq / (eta(i,k) + etah * c0 * dz)
847 qcko(i,k) = qlk + qrch
848 pwo(i,k) = etah * c0 * dz * qlk
849 cnvwt(i,k) = etah * qlk * g / dp
856 c exchange ktcon with ktcon1
866 c this section is ready for cloud water
870 c compute liquid and vapor separation at cloud top
876 gamma = el2orc * qeso(i,k) / (to(i,k)**2)
878 & + gamma * dbyo(i,k) / (hvap * (1. + gamma))
879 dq = qcko(i,k) - qrch
881 c check
if there is excess moisture to release latent heat
891 c--- compute precipitation efficiency in terms of windshear
906 if(k.gt.kb(i).and.k.le.ktcon(i))
then
907 shear= sqrt((uo(i,k)-uo(i,k-1)) ** 2
908 & + (vo(i,k)-vo(i,k-1)) ** 2)
909 vshear(i) = vshear(i) + shear
916 vshear(i) = 1.e3 * vshear(i) / (zi(i,ktcon(i))-zi(i,kb(i)))
917 e1=1.591-.639*vshear(i)
918 & +.0953*(vshear(i)**2)-.00496*(vshear(i)**3)
921 edt(i) = min(edt(i),val)
923 edt(i) = max(edt(i),val)
927 c--- what would the change be, that a cloud with unit mass
928 c--- will
do to the environment?
934 if(cnvflg(i) .and. k .le. kmax(i))
then
943 c--- changed due to subsidence and entrainment
948 if(k.gt.kb(i).and.k.lt.ktcon(i))
then
949 dp = 1000. * del(i,k)
950 dz = zi(i,k) - zi(i,k-1)
953 dv2h = .5 * (heo(i,k) + heo(i,k-1))
956 dv2q = .5 * (qo(i,k) + qo(i,k-1))
959 dv2u = .5 * (uo(i,k) + uo(i,k-1))
962 dv2v = .5 * (vo(i,k) + vo(i,k-1))
965 tem = 0.5 * (xlamue(i,k)+xlamue(i,k-1))
968 dellah(i,k) = dellah(i,k) +
969 & ( eta(i,k)*dv1h - eta(i,k-1)*dv3h
970 & - tem*eta(i,k-1)*dv2h*dz
971 & + tem1*eta(i,k-1)*.5*(hcko(i,k)+hcko(i,k-1))*dz
974 dellaq(i,k) = dellaq(i,k) +
975 & ( eta(i,k)*dv1q - eta(i,k-1)*dv3q
976 & - tem*eta(i,k-1)*dv2q*dz
977 & + tem1*eta(i,k-1)*.5*(qrcko(i,k)+qcko(i,k-1))*dz
980 dellau(i,k) = dellau(i,k) +
981 & ( eta(i,k)*dv1u - eta(i,k-1)*dv3u
982 & - tem*eta(i,k-1)*dv2u*dz
983 & + tem1*eta(i,k-1)*.5*(ucko(i,k)+ucko(i,k-1))*dz
984 & - pgcon*eta(i,k-1)*(dv1u-dv3u)
987 dellav(i,k) = dellav(i,k) +
988 & ( eta(i,k)*dv1v - eta(i,k-1)*dv3v
989 & - tem*eta(i,k-1)*dv2v*dz
990 & + tem1*eta(i,k-1)*.5*(vcko(i,k)+vcko(i,k-1))*dz
991 & - pgcon*eta(i,k-1)*(dv1v-dv3v)
1004 dp = 1000. * del(i,indx)
1005 dv1h = heo(i,indx-1)
1006 dellah(i,indx) = eta(i,indx-1) *
1007 & (hcko(i,indx-1) - dv1h) * g / dp
1009 dellaq(i,indx) = eta(i,indx-1) *
1010 & (qcko(i,indx-1) - dv1q) * g / dp
1012 dellau(i,indx) = eta(i,indx-1) *
1013 & (ucko(i,indx-1) - dv1u) * g / dp
1015 dellav(i,indx) = eta(i,indx-1) *
1016 & (vcko(i,indx-1) - dv1v) * g / dp
1020 dellal(i,indx) = eta(i,indx-1) *
1021 & qlko_ktcon(i) * g / dp
1025 c mass flux at cloud base for shallow convection
1037 ptem = g*sflx(i)*hpbl(i)/t1(i,1)
1039 tem = po(i,k)*100. / (rd*t1(i,k))
1040 xmb(i) = betaw*tem*wstar(i)
1041 xmb(i) = min(xmb(i),xmbmax(i))
1051 if (cnvflg(i) .and. k .le. kmax(i))
then
1052 qeso(i,k) = 0.01 * fpvs(t1(i,k))
1053 qeso(i,k) = eps * qeso(i,k) / (pfld(i,k) + epsm1*qeso(i,k))
1055 qeso(i,k) = max(qeso(i,k), val )
1075 if(k.gt.kb(i).and.k.le.ktcon(i))
then
1076 dellat = (dellah(i,k) - hvap * dellaq(i,k)) / cp
1077 t1(i,k) = t1(i,k) + dellat * xmb(i) * dt2
1078 q1(i,k) = q1(i,k) + dellaq(i,k) * xmb(i) * dt2
1082 u1(i,k) = u1(i,k) + dellau(i,k) * xmb(i) * dt2
1083 v1(i,k) = v1(i,k) + dellav(i,k) * xmb(i) * dt2
1084 dp = 1000. * del(i,k)
1085 delhbar(i) = delhbar(i) + dellah(i,k)*xmb(i)*dp/g
1086 delqbar(i) = delqbar(i) + dellaq(i,k)*xmb(i)*dp/g
1087 deltbar(i) = deltbar(i) + dellat*xmb(i)*dp/g
1088 delubar(i) = delubar(i) + dellau(i,k)*xmb(i)*dp/g
1089 delvbar(i) = delvbar(i) + dellav(i,k)*xmb(i)*dp/g
1098 if(k.gt.kb(i).and.k.le.ktcon(i))
then
1099 qeso(i,k) = 0.01 * fpvs(t1(i,k))
1100 qeso(i,k) = eps * qeso(i,k)/(pfld(i,k) + epsm1*qeso(i,k))
1102 qeso(i,k) = max(qeso(i,k), val )
1118 if(k.lt.ktcon(i).and.k.gt.kb(i))
then
1119 rntot(i) = rntot(i) + pwo(i,k) * xmb(i) * .001 * dt2
1132 if (k .le. kmax(i))
then
1137 if(k.lt.ktcon(i).and.k.gt.kb(i))
then
1138 rn(i) = rn(i) + pwo(i,k) * xmb(i) * .001 * dt2
1141 if(flg(i).and.k.lt.ktcon(i))
then
1142 evef = edt(i) * evfact
1143 if(islimsk(i) == 1) evef=edt(i) * evfactl
1145 c
if(islimsk(i) == 1) evef = 0.
1146 qcond(i) = evef * (q1(i,k) - qeso(i,k))
1147 & / (1. + el2orc * qeso(i,k) / t1(i,k)**2)
1148 dp = 1000. * del(i,k)
1149 if(rn(i).gt.0..and.qcond(i).lt.0.)
then
1150 qevap(i) = -qcond(i) * (1.-exp(-.32*sqrt(dt2*rn(i))))
1151 qevap(i) = min(qevap(i), rn(i)*1000.*g/dp)
1152 delq2(i) = delqev(i) + .001 * qevap(i) * dp / g
1154 if(rn(i).gt.0..and.qcond(i).lt.0..and.
1155 & delq2(i).gt.rntot(i))
then
1156 qevap(i) = 1000.* g * (rntot(i) - delqev(i)) / dp
1159 if(rn(i).gt.0..and.qevap(i).gt.0.)
then
1161 tem1 = qevap(i) * tem
1162 if(tem1.gt.rn(i))
then
1163 qevap(i) = rn(i) / tem
1166 rn(i) = rn(i) - tem1
1168 q1(i,k) = q1(i,k) + qevap(i)
1169 t1(i,k) = t1(i,k) - elocp * qevap(i)
1170 deltv(i) = - elocp*qevap(i)/dt2
1171 delq(i) = + qevap(i)/dt2
1172 delqev(i) = delqev(i) + .001*dp*qevap(i)/g
1174 dellaq(i,k) = dellaq(i,k) + delq(i) / xmb(i)
1175 delqbar(i) = delqbar(i) + delq(i)*dp/g
1176 deltbar(i) = deltbar(i) + deltv(i)*dp/g
1195 if(rn(i).lt.0..or..not.flg(i)) rn(i) = 0.
1202 c convective cloud water
1207 if (cnvflg(i) .and. rn(i).gt.0.)
then
1208 if (k.ge.kbcon(i).and.k.lt.ktcon(i))
then
1209 cnvw(i,k) = cnvwt(i,k) * xmb(i) * dt2
1216 c convective cloud cover
1221 if (cnvflg(i) .and. rn(i).gt.0.)
then
1222 if (k.ge.kbcon(i).and.k.lt.ktcon(i))
then
1223 cnvc(i,k) = 0.04 * log(1. + 675. * eta(i,k) * xmb(i))
1224 cnvc(i,k) = min(cnvc(i,k), 0.2)
1225 cnvc(i,k) = max(cnvc(i,k), 0.0)
1235 if (ncloud.gt.0)
then
1240 if (k.gt.kb(i).and.k.le.ktcon(i))
then
1241 tem = dellal(i,k) * xmb(i) * dt2
1242 tem1 = max(0.0, min(1.0, (tcr-t1(i,k))*tcrf))
1243 if (ql(i,k,2) .gt. -999.0)
then
1244 ql(i,k,1) = ql(i,k,1) + tem * tem1
1245 ql(i,k,2) = ql(i,k,2) + tem *(1.0-tem1)
1247 ql(i,k,1) = ql(i,k,1) + tem
1262 if(k.ge.kb(i) .and. k.lt.ktop(i))
then
1263 ud_mf(i,k) = eta(i,k) * xmb(i) * dt2
1272 dt_mf(i,k) = ud_mf(i,k)
real(kind=kind_phys), parameter con_cliq
spec heat H2O liq ( )
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 ( )
subroutine shalcnv(im, ix, km, jcap, delt, delp, prslp, psp, phil, ql, q1, t1, u1, v1, rn, kbot, ktop, kcnv, islimsk, dot, ncloud, hpbl, heat, evap, ud_mf, dt_mf, cnvw, cnvc)
This subroutine contains the entirety of the SAS shallow convection scheme.
real(kind=kind_phys), parameter con_cp
spec heat air at p ( )
real(kind=kind_phys), parameter con_hvap
lat heat H2O cond ( )
real(kind=kind_phys), parameter con_rd
gas constant air ( )