1 |
module cv3_unsat_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
SUBROUTINE cv3_unsat(nloc,ncum,nd,na,icb,inb & |
8 |
,t,rr,rs,gz,u,v,p,ph & |
9 |
,th,tv,lv,cpn,ep,sigp,clw & |
10 |
,m,ment,elij,delt,plcl & |
11 |
,mp,rp,up,vp,wt,water,evap,b) |
12 |
use cv3_param_m |
13 |
use cvthermo |
14 |
use cvflag |
15 |
|
16 |
|
17 |
! inputs: |
18 |
integer, intent(in):: ncum, nd, na, nloc |
19 |
integer icb(nloc), inb(nloc) |
20 |
real, intent(in):: delt |
21 |
real plcl(nloc) |
22 |
real t(nloc,nd), rr(nloc,nd), rs(nloc,nd) |
23 |
real u(nloc,nd), v(nloc,nd) |
24 |
real p(nloc,nd), ph(nloc,nd+1) |
25 |
real th(nloc,na), gz(nloc,na) |
26 |
real lv(nloc,na), ep(nloc,na), sigp(nloc,na), clw(nloc,na) |
27 |
real cpn(nloc,na), tv(nloc,na) |
28 |
real m(nloc,na), ment(nloc,na,na), elij(nloc,na,na) |
29 |
|
30 |
! outputs: |
31 |
real mp(nloc,na), rp(nloc,na), up(nloc,na), vp(nloc,na) |
32 |
real water(nloc,na), evap(nloc,na), wt(nloc,na) |
33 |
real b(nloc,na) |
34 |
|
35 |
! local variables |
36 |
integer i,j,il,num1 |
37 |
real tinv, delti |
38 |
real awat, afac, afac1, afac2, bfac |
39 |
real pr1, pr2, sigt, b6, c6, revap, tevap, delth |
40 |
real amfac, amp2, xf, tf, fac2, ur, sru, fac, d, af, bf |
41 |
real ampmax |
42 |
real lvcp(nloc,na) |
43 |
real wdtrain(nloc) |
44 |
logical lwork(nloc) |
45 |
|
46 |
|
47 |
!------------------------------------------------------ |
48 |
|
49 |
delti = 1./delt |
50 |
tinv=1./3. |
51 |
|
52 |
mp(:,:)=0. |
53 |
|
54 |
do i=1,nl |
55 |
do il=1,ncum |
56 |
mp(il,i)=0.0 |
57 |
rp(il,i)=rr(il,i) |
58 |
up(il,i)=u(il,i) |
59 |
vp(il,i)=v(il,i) |
60 |
wt(il,i)=0.001 |
61 |
water(il,i)=0.0 |
62 |
evap(il,i)=0.0 |
63 |
b(il,i)=0.0 |
64 |
lvcp(il,i)=lv(il,i)/cpn(il,i) |
65 |
enddo |
66 |
enddo |
67 |
|
68 |
! |
69 |
! *** check whether ep(inb)=0, if so, skip precipitating *** |
70 |
! *** downdraft calculation *** |
71 |
! |
72 |
|
73 |
do il=1,ncum |
74 |
lwork(il)=.TRUE. |
75 |
if(ep(il,inb(il)).lt.0.0001)lwork(il)=.FALSE. |
76 |
enddo |
77 |
|
78 |
call zilch(wdtrain,ncum) |
79 |
|
80 |
DO i=nl+1,1,-1 |
81 |
|
82 |
num1=0 |
83 |
do il=1,ncum |
84 |
if ( i.le.inb(il) .and. lwork(il) ) num1=num1+1 |
85 |
enddo |
86 |
if (num1.le.0) cycle |
87 |
|
88 |
! |
89 |
! *** integrate liquid water equation to find condensed water *** |
90 |
! *** and condensed water flux *** |
91 |
! |
92 |
|
93 |
! |
94 |
! *** begin downdraft loop *** |
95 |
! |
96 |
|
97 |
! |
98 |
! *** calculate detrained precipitation *** |
99 |
! |
100 |
do il=1,ncum |
101 |
if (i.le.inb(il) .and. lwork(il)) then |
102 |
if (cvflag_grav) then |
103 |
wdtrain(il)=grav*ep(il,i)*m(il,i)*clw(il,i) |
104 |
else |
105 |
wdtrain(il)=10.0*ep(il,i)*m(il,i)*clw(il,i) |
106 |
endif |
107 |
endif |
108 |
enddo |
109 |
|
110 |
if(i.gt.1)then |
111 |
do j=1,i-1 |
112 |
do il=1,ncum |
113 |
if (i.le.inb(il) .and. lwork(il)) then |
114 |
awat=elij(il,j,i)-(1.-ep(il,i))*clw(il,i) |
115 |
awat=amax1(awat,0.0) |
116 |
if (cvflag_grav) then |
117 |
wdtrain(il)=wdtrain(il)+grav*awat*ment(il,j,i) |
118 |
else |
119 |
wdtrain(il)=wdtrain(il)+10.0*awat*ment(il,j,i) |
120 |
endif |
121 |
endif |
122 |
enddo |
123 |
end do |
124 |
endif |
125 |
|
126 |
! |
127 |
! *** find rain water and evaporation using provisional *** |
128 |
! *** estimates of rp(i)and rp(i-1) *** |
129 |
! |
130 |
|
131 |
do il=1,ncum |
132 |
|
133 |
if (i.le.inb(il) .and. lwork(il)) then |
134 |
|
135 |
wt(il,i)=45.0 |
136 |
|
137 |
if(i.lt.inb(il))then |
138 |
rp(il,i)=rp(il,i+1) & |
139 |
+(cpd*(t(il,i+1)-t(il,i))+gz(il,i+1)-gz(il,i))/lv(il,i) |
140 |
rp(il,i)=0.5*(rp(il,i)+rr(il,i)) |
141 |
endif |
142 |
rp(il,i)=amax1(rp(il,i),0.0) |
143 |
rp(il,i)=amin1(rp(il,i),rs(il,i)) |
144 |
rp(il,inb(il))=rr(il,inb(il)) |
145 |
|
146 |
if(i.eq.1)then |
147 |
afac=p(il,1)*(rs(il,1)-rp(il,1))/(1.0e4+2000.0*p(il,1)*rs(il,1)) |
148 |
else |
149 |
rp(il,i-1)=rp(il,i) & |
150 |
+(cpd*(t(il,i)-t(il,i-1))+gz(il,i)-gz(il,i-1))/lv(il,i) |
151 |
rp(il,i-1)=0.5*(rp(il,i-1)+rr(il,i-1)) |
152 |
rp(il,i-1)=amin1(rp(il,i-1),rs(il,i-1)) |
153 |
rp(il,i-1)=amax1(rp(il,i-1),0.0) |
154 |
afac1=p(il,i)*(rs(il,i)-rp(il,i))/(1.0e4+2000.0*p(il,i)*rs(il,i)) |
155 |
afac2=p(il,i-1)*(rs(il,i-1)-rp(il,i-1)) & |
156 |
/(1.0e4+2000.0*p(il,i-1)*rs(il,i-1)) |
157 |
afac=0.5*(afac1+afac2) |
158 |
endif |
159 |
if(i.eq.inb(il))afac=0.0 |
160 |
afac=amax1(afac,0.0) |
161 |
bfac=1./(sigd*wt(il,i)) |
162 |
! |
163 |
!jyg1 |
164 |
!cc sigt=1.0 |
165 |
!cc if(i.ge.icb)sigt=sigp(i) |
166 |
! prise en compte de la variation progressive de sigt dans |
167 |
! les couches icb et icb-1: |
168 |
! pour plcl<ph(i+1), pr1=0 & pr2=1 |
169 |
! pour plcl>ph(i), pr1=1 & pr2=0 |
170 |
! pour ph(i+1)<plcl<ph(i), pr1 est la proportion a cheval |
171 |
! sur le nuage, et pr2 est la proportion sous la base du |
172 |
! nuage. |
173 |
pr1=(plcl(il)-ph(il,i+1))/(ph(il,i)-ph(il,i+1)) |
174 |
pr1=max(0.,min(1.,pr1)) |
175 |
pr2=(ph(il,i)-plcl(il))/(ph(il,i)-ph(il,i+1)) |
176 |
pr2=max(0.,min(1.,pr2)) |
177 |
sigt=sigp(il,i)*pr1+pr2 |
178 |
!jyg2 |
179 |
! |
180 |
b6=bfac*50.*sigd*(ph(il,i)-ph(il,i+1))*sigt*afac |
181 |
c6=water(il,i+1)+bfac*wdtrain(il) & |
182 |
-50.*sigd*bfac*(ph(il,i)-ph(il,i+1))*evap(il,i+1) |
183 |
if(c6.gt.0.0)then |
184 |
revap=0.5*(-b6+sqrt(b6*b6+4.*c6)) |
185 |
evap(il,i)=sigt*afac*revap |
186 |
water(il,i)=revap*revap |
187 |
else |
188 |
evap(il,i)=-evap(il,i+1) & |
189 |
+0.02*(wdtrain(il)+sigd*wt(il,i)*water(il,i+1)) & |
190 |
/(sigd*(ph(il,i)-ph(il,i+1))) |
191 |
end if |
192 |
! |
193 |
! *** calculate precipitating downdraft mass flux under *** |
194 |
! *** hydrostatic approximation *** |
195 |
! |
196 |
if (i.ne.1) then |
197 |
|
198 |
tevap=amax1(0.0,evap(il,i)) |
199 |
delth=amax1(0.001,(th(il,i)-th(il,i-1))) |
200 |
if (cvflag_grav) then |
201 |
mp(il,i)=100.*ginv*lvcp(il,i)*sigd*tevap & |
202 |
*(p(il,i-1)-p(il,i))/delth |
203 |
else |
204 |
mp(il,i)=10.*lvcp(il,i)*sigd*tevap*(p(il,i-1)-p(il,i))/delth |
205 |
endif |
206 |
! |
207 |
! *** if hydrostatic assumption fails, *** |
208 |
! *** solve cubic difference equation for downdraft theta *** |
209 |
! *** and mass flux from two simultaneous differential eqns *** |
210 |
! |
211 |
amfac=sigd*sigd*70.0*ph(il,i)*(p(il,i-1)-p(il,i)) & |
212 |
*(th(il,i)-th(il,i-1))/(tv(il,i)*th(il,i)) |
213 |
amp2=abs(mp(il,i+1)*mp(il,i+1)-mp(il,i)*mp(il,i)) |
214 |
if(amp2.gt.(0.1*amfac))then |
215 |
xf=100.0*sigd*sigd*sigd*(ph(il,i)-ph(il,i+1)) |
216 |
tf=b(il,i)-5.0*(th(il,i)-th(il,i-1))*t(il,i) & |
217 |
/(lvcp(il,i)*sigd*th(il,i)) |
218 |
af=xf*tf+mp(il,i+1)*mp(il,i+1)*tinv |
219 |
bf=2.*(tinv*mp(il,i+1))**3+tinv*mp(il,i+1)*xf*tf & |
220 |
+50.*(p(il,i-1)-p(il,i))*xf*tevap |
221 |
fac2=1.0 |
222 |
if(bf.lt.0.0)fac2=-1.0 |
223 |
bf=abs(bf) |
224 |
ur=0.25*bf*bf-af*af*af*tinv*tinv*tinv |
225 |
if(ur.ge.0.0)then |
226 |
sru=sqrt(ur) |
227 |
fac=1.0 |
228 |
if((0.5*bf-sru).lt.0.0)fac=-1.0 |
229 |
mp(il,i)=mp(il,i+1)*tinv+(0.5*bf+sru)**tinv & |
230 |
+fac*(abs(0.5*bf-sru))**tinv |
231 |
else |
232 |
d=atan(2.*sqrt(-ur)/(bf+1.0e-28)) |
233 |
if(fac2.lt.0.0)d=3.14159-d |
234 |
mp(il,i)=mp(il,i+1)*tinv+2.*sqrt(af*tinv)*cos(d*tinv) |
235 |
endif |
236 |
mp(il,i)=amax1(0.0,mp(il,i)) |
237 |
|
238 |
if (cvflag_grav) then |
239 |
!jyg : il y a vraisemblablement une erreur dans la ligne 2 suivante: |
240 |
! il faut diviser par (mp(il,i)*sigd*grav) et non par (mp(il,i)+sigd*0.1). |
241 |
! Et il faut bien revoir les facteurs 100. |
242 |
b(il,i-1)=b(il,i)+100.0*(p(il,i-1)-p(il,i))*tevap & |
243 |
/(mp(il,i)+sigd*0.1) & |
244 |
-10.0*(th(il,i)-th(il,i-1))*t(il,i)/(lvcp(il,i)*sigd*th(il,i)) |
245 |
else |
246 |
b(il,i-1)=b(il,i)+100.0*(p(il,i-1)-p(il,i))*tevap & |
247 |
/(mp(il,i)+sigd*0.1) & |
248 |
-10.0*(th(il,i)-th(il,i-1))*t(il,i)/(lvcp(il,i)*sigd*th(il,i)) |
249 |
endif |
250 |
b(il,i-1)=amax1(b(il,i-1),0.0) |
251 |
endif |
252 |
! |
253 |
! *** limit magnitude of mp(i) to meet cfl condition *** |
254 |
! |
255 |
ampmax=2.0*(ph(il,i)-ph(il,i+1))*delti |
256 |
amp2=2.0*(ph(il,i-1)-ph(il,i))*delti |
257 |
ampmax=amin1(ampmax,amp2) |
258 |
mp(il,i)=amin1(mp(il,i),ampmax) |
259 |
! |
260 |
! *** force mp to decrease linearly to zero *** |
261 |
! *** between cloud base and the surface *** |
262 |
! |
263 |
if(p(il,i).gt.p(il,icb(il)))then |
264 |
mp(il,i)=mp(il,icb(il))*(p(il,1)-p(il,i))/(p(il,1)-p(il,icb(il))) |
265 |
endif |
266 |
|
267 |
endif ! i.eq.1 |
268 |
! |
269 |
! *** find mixing ratio of precipitating downdraft *** |
270 |
! |
271 |
|
272 |
if (i.ne.inb(il)) then |
273 |
|
274 |
rp(il,i)=rr(il,i) |
275 |
|
276 |
if(mp(il,i).gt.mp(il,i+1))then |
277 |
|
278 |
if (cvflag_grav) then |
279 |
rp(il,i)=rp(il,i+1)*mp(il,i+1)+rr(il,i)*(mp(il,i)-mp(il,i+1)) & |
280 |
+100.*ginv*0.5*sigd*(ph(il,i)-ph(il,i+1)) & |
281 |
*(evap(il,i+1)+evap(il,i)) |
282 |
else |
283 |
rp(il,i)=rp(il,i+1)*mp(il,i+1)+rr(il,i)*(mp(il,i)-mp(il,i+1)) & |
284 |
+5.*sigd*(ph(il,i)-ph(il,i+1)) & |
285 |
*(evap(il,i+1)+evap(il,i)) |
286 |
endif |
287 |
rp(il,i)=rp(il,i)/mp(il,i) |
288 |
up(il,i)=up(il,i+1)*mp(il,i+1)+u(il,i)*(mp(il,i)-mp(il,i+1)) |
289 |
up(il,i)=up(il,i)/mp(il,i) |
290 |
vp(il,i)=vp(il,i+1)*mp(il,i+1)+v(il,i)*(mp(il,i)-mp(il,i+1)) |
291 |
vp(il,i)=vp(il,i)/mp(il,i) |
292 |
|
293 |
else |
294 |
|
295 |
if(mp(il,i+1).gt.1.0e-16)then |
296 |
if (cvflag_grav) then |
297 |
rp(il,i)=rp(il,i+1) & |
298 |
+100.*ginv*0.5*sigd*(ph(il,i)-ph(il,i+1)) & |
299 |
*(evap(il,i+1)+evap(il,i))/mp(il,i+1) |
300 |
else |
301 |
rp(il,i)=rp(il,i+1) & |
302 |
+5.*sigd*(ph(il,i)-ph(il,i+1)) & |
303 |
*(evap(il,i+1)+evap(il,i))/mp(il,i+1) |
304 |
endif |
305 |
up(il,i)=up(il,i+1) |
306 |
vp(il,i)=vp(il,i+1) |
307 |
|
308 |
endif |
309 |
endif |
310 |
rp(il,i)=amin1(rp(il,i),rs(il,i)) |
311 |
rp(il,i)=amax1(rp(il,i),0.0) |
312 |
|
313 |
endif |
314 |
endif |
315 |
end do |
316 |
|
317 |
end DO |
318 |
|
319 |
end SUBROUTINE cv3_unsat |
320 |
|
321 |
end module cv3_unsat_m |