1 |
|
2 |
SUBROUTINE cv3_undilute2(nloc,ncum,nd,icb,icbs,nk & |
3 |
,tnk,qnk,gznk,t,q,qs,gz & |
4 |
,p,h,tv,lv,pbase,buoybase,plcl & |
5 |
,inb,tp,tvp,clw,hp,ep,sigp,buoy) |
6 |
use conema3_m |
7 |
use cv3_param_m |
8 |
use cvthermo |
9 |
implicit none |
10 |
|
11 |
!--------------------------------------------------------------------- |
12 |
! Purpose: |
13 |
! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
14 |
! & |
15 |
! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE |
16 |
! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD |
17 |
! & |
18 |
! FIND THE LEVEL OF NEUTRAL BUOYANCY |
19 |
! |
20 |
! Main differences convect3/convect4: |
21 |
! - icbs (input) is the first level above LCL (may differ from icb) |
22 |
! - many minor differences in the iterations |
23 |
! - condensed water not removed from tvp in convect3 |
24 |
! - vertical profile of buoyancy computed here (use of buoybase) |
25 |
! - the determination of inb is different |
26 |
! - no inb1, only inb in output |
27 |
!--------------------------------------------------------------------- |
28 |
|
29 |
|
30 |
! inputs: |
31 |
integer ncum, nd, nloc |
32 |
integer icb(nloc), icbs(nloc), nk(nloc) |
33 |
real t(nloc,nd), q(nloc,nd), qs(nloc,nd), gz(nloc,nd) |
34 |
real p(nloc,nd) |
35 |
real tnk(nloc), qnk(nloc), gznk(nloc) |
36 |
real lv(nloc,nd), tv(nloc,nd), h(nloc,nd) |
37 |
real pbase(nloc), buoybase(nloc), plcl(nloc) |
38 |
|
39 |
! outputs: |
40 |
integer inb(nloc) |
41 |
real tp(nloc,nd), tvp(nloc,nd), clw(nloc,nd) |
42 |
real ep(nloc,nd), sigp(nloc,nd), hp(nloc,nd) |
43 |
real buoy(nloc,nd) |
44 |
|
45 |
! local variables: |
46 |
integer i, k |
47 |
real tg,qg,ahg,alv,s,tc,es,denom,rg,tca,elacrit |
48 |
real by, defrac, pden |
49 |
real ah0(nloc), cape(nloc), capem(nloc), byp(nloc) |
50 |
logical lcape(nloc) |
51 |
|
52 |
!===================================================================== |
53 |
! --- SOME INITIALIZATIONS |
54 |
!===================================================================== |
55 |
|
56 |
do 170 k=1,nl |
57 |
do 160 i=1,ncum |
58 |
ep(i,k)=0.0 |
59 |
sigp(i,k)=spfac |
60 |
160 continue |
61 |
170 continue |
62 |
|
63 |
!===================================================================== |
64 |
! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
65 |
!===================================================================== |
66 |
! |
67 |
! --- The procedure is to solve the equation. |
68 |
! cp*tp+L*qp+phi=cp*tnk+L*qnk+gznk. |
69 |
! |
70 |
! *** Calculate certain parcel quantities, including static energy *** |
71 |
! |
72 |
! |
73 |
do 240 i=1,ncum |
74 |
ah0(i)=(cpd*(1.-qnk(i))+cl*qnk(i))*tnk(i) & |
75 |
+qnk(i)*(lv0-clmcpv*(tnk(i)-273.15))+gznk(i) |
76 |
240 continue |
77 |
! |
78 |
! |
79 |
! *** Find lifted parcel quantities above cloud base *** |
80 |
! |
81 |
! |
82 |
do 300 k=minorig+1,nl |
83 |
do 290 i=1,ncum |
84 |
! ori if(k.ge.(icb(i)+1))then |
85 |
if(k.ge.(icbs(i)+1))then ! convect3 |
86 |
tg=t(i,k) |
87 |
qg=qs(i,k) |
88 |
!debug alv=lv0-clmcpv*(t(i,k)-t0) |
89 |
alv=lv0-clmcpv*(t(i,k)-273.15) |
90 |
! |
91 |
! First iteration. |
92 |
! |
93 |
! ori s=cpd+alv*alv*qg/(rrv*t(i,k)*t(i,k)) |
94 |
s=cpd*(1.-qnk(i))+cl*qnk(i) & |
95 |
+alv*alv*qg/(rrv*t(i,k)*t(i,k)) ! convect3 |
96 |
s=1./s |
97 |
! ori ahg=cpd*tg+(cl-cpd)*qnk(i)*t(i,k)+alv*qg+gz(i,k) |
98 |
ahg=cpd*tg+(cl-cpd)*qnk(i)*tg+alv*qg+gz(i,k) ! convect3 |
99 |
tg=tg+s*(ah0(i)-ahg) |
100 |
! ori tg=max(tg,35.0) |
101 |
!debug tc=tg-t0 |
102 |
tc=tg-273.15 |
103 |
denom=243.5+tc |
104 |
denom=MAX(denom,1.0) ! convect3 |
105 |
! ori if(tc.ge.0.0)then |
106 |
es=6.112*exp(17.67*tc/denom) |
107 |
! ori else |
108 |
! ori es=exp(23.33086-6111.72784/tg+0.15215*log(tg)) |
109 |
! ori endif |
110 |
qg=eps*es/(p(i,k)-es*(1.-eps)) |
111 |
! |
112 |
! Second iteration. |
113 |
! |
114 |
! ori s=cpd+alv*alv*qg/(rrv*t(i,k)*t(i,k)) |
115 |
! ori s=1./s |
116 |
! ori ahg=cpd*tg+(cl-cpd)*qnk(i)*t(i,k)+alv*qg+gz(i,k) |
117 |
ahg=cpd*tg+(cl-cpd)*qnk(i)*tg+alv*qg+gz(i,k) ! convect3 |
118 |
tg=tg+s*(ah0(i)-ahg) |
119 |
! ori tg=max(tg,35.0) |
120 |
!debug tc=tg-t0 |
121 |
tc=tg-273.15 |
122 |
denom=243.5+tc |
123 |
denom=MAX(denom,1.0) ! convect3 |
124 |
! ori if(tc.ge.0.0)then |
125 |
es=6.112*exp(17.67*tc/denom) |
126 |
! ori else |
127 |
! ori es=exp(23.33086-6111.72784/tg+0.15215*log(tg)) |
128 |
! ori endif |
129 |
qg=eps*es/(p(i,k)-es*(1.-eps)) |
130 |
! |
131 |
!debug alv=lv0-clmcpv*(t(i,k)-t0) |
132 |
alv=lv0-clmcpv*(t(i,k)-273.15) |
133 |
! print*,'cpd dans convect2 ',cpd |
134 |
! print*,'tp(i,k),ah0(i),cl,cpd,qnk(i),t(i,k),gz(i,k),alv,qg,cpd' |
135 |
! print*,tp(i,k),ah0(i),cl,cpd,qnk(i),t(i,k),gz(i,k),alv,qg,cpd |
136 |
|
137 |
! ori c approximation here: |
138 |
! ori tp(i,k)=(ah0(i)-(cl-cpd)*qnk(i)*t(i,k)-gz(i,k)-alv*qg)/cpd |
139 |
|
140 |
! convect3: no approximation: |
141 |
tp(i,k)=(ah0(i)-gz(i,k)-alv*qg)/(cpd+(cl-cpd)*qnk(i)) |
142 |
|
143 |
clw(i,k)=qnk(i)-qg |
144 |
clw(i,k)=max(0.0,clw(i,k)) |
145 |
rg=qg/(1.-qnk(i)) |
146 |
! ori tvp(i,k)=tp(i,k)*(1.+rg*epsi) |
147 |
! convect3: (qg utilise au lieu du vrai mixing ratio rg): |
148 |
tvp(i,k)=tp(i,k)*(1.+qg/eps-qnk(i)) ! whole thing |
149 |
endif |
150 |
290 continue |
151 |
300 continue |
152 |
! |
153 |
!===================================================================== |
154 |
! --- SET THE PRECIPITATION EFFICIENCIES AND THE FRACTION OF |
155 |
! --- PRECIPITATION FALLING OUTSIDE OF CLOUD |
156 |
! --- THESE MAY BE FUNCTIONS OF TP(I), P(I) AND CLW(I) |
157 |
!===================================================================== |
158 |
! |
159 |
! ori do 320 k=minorig+1,nl |
160 |
do 320 k=1,nl ! convect3 |
161 |
do 310 i=1,ncum |
162 |
pden=ptcrit-pbcrit |
163 |
ep(i,k)=(plcl(i)-p(i,k)-pbcrit)/pden*epmax |
164 |
ep(i,k)=amax1(ep(i,k),0.0) |
165 |
ep(i,k)=amin1(ep(i,k),epmax) |
166 |
sigp(i,k)=spfac |
167 |
! ori if(k.ge.(nk(i)+1))then |
168 |
! ori tca=tp(i,k)-t0 |
169 |
! ori if(tca.ge.0.0)then |
170 |
! ori elacrit=elcrit |
171 |
! ori else |
172 |
! ori elacrit=elcrit*(1.0-tca/tlcrit) |
173 |
! ori endif |
174 |
! ori elacrit=max(elacrit,0.0) |
175 |
! ori ep(i,k)=1.0-elacrit/max(clw(i,k),1.0e-8) |
176 |
! ori ep(i,k)=max(ep(i,k),0.0 ) |
177 |
! ori ep(i,k)=min(ep(i,k),1.0 ) |
178 |
! ori sigp(i,k)=sigs |
179 |
! ori endif |
180 |
310 continue |
181 |
320 continue |
182 |
! |
183 |
!===================================================================== |
184 |
! --- CALCULATE VIRTUAL TEMPERATURE AND LIFTED PARCEL |
185 |
! --- VIRTUAL TEMPERATURE |
186 |
!===================================================================== |
187 |
! |
188 |
! dans convect3, tvp est calcule en une seule fois, et sans retirer |
189 |
! l'eau condensee (~> reversible CAPE) |
190 |
! |
191 |
! ori do 340 k=minorig+1,nl |
192 |
! ori do 330 i=1,ncum |
193 |
! ori if(k.ge.(icb(i)+1))then |
194 |
! ori tvp(i,k)=tvp(i,k)*(1.0-qnk(i)+ep(i,k)*clw(i,k)) |
195 |
! oric print*,'i,k,tvp(i,k),qnk(i),ep(i,k),clw(i,k)' |
196 |
! oric print*, i,k,tvp(i,k),qnk(i),ep(i,k),clw(i,k) |
197 |
! ori endif |
198 |
! ori 330 continue |
199 |
! ori 340 continue |
200 |
|
201 |
! ori do 350 i=1,ncum |
202 |
! ori tvp(i,nlp)=tvp(i,nl)-(gz(i,nlp)-gz(i,nl))/cpd |
203 |
! ori 350 continue |
204 |
|
205 |
do 350 i=1,ncum ! convect3 |
206 |
tp(i,nlp)=tp(i,nl) ! convect3 |
207 |
350 continue ! convect3 |
208 |
! |
209 |
!===================================================================== |
210 |
! --- EFFECTIVE VERTICAL PROFILE OF BUOYANCY (convect3 only): |
211 |
!===================================================================== |
212 |
|
213 |
!-- this is for convect3 only: |
214 |
|
215 |
! first estimate of buoyancy: |
216 |
|
217 |
do 500 i=1,ncum |
218 |
do 501 k=1,nl |
219 |
buoy(i,k)=tvp(i,k)-tv(i,k) |
220 |
501 continue |
221 |
500 continue |
222 |
|
223 |
! set buoyancy=buoybase for all levels below base |
224 |
! for safety, set buoy(icb)=buoybase |
225 |
|
226 |
do 505 i=1,ncum |
227 |
do 506 k=1,nl |
228 |
if((k.ge.icb(i)).and.(k.le.nl).and.(p(i,k).ge.pbase(i)))then |
229 |
buoy(i,k)=buoybase(i) |
230 |
endif |
231 |
506 continue |
232 |
buoy(icb(i),k)=buoybase(i) |
233 |
505 continue |
234 |
|
235 |
!-- end convect3 |
236 |
|
237 |
!===================================================================== |
238 |
! --- FIND THE FIRST MODEL LEVEL (INB) ABOVE THE PARCEL'S |
239 |
! --- LEVEL OF NEUTRAL BUOYANCY |
240 |
!===================================================================== |
241 |
! |
242 |
!-- this is for convect3 only: |
243 |
|
244 |
do 510 i=1,ncum |
245 |
inb(i)=nl-1 |
246 |
510 continue |
247 |
|
248 |
do 530 i=1,ncum |
249 |
do 535 k=1,nl-1 |
250 |
if ((k.ge.icb(i)).and.(buoy(i,k).lt.dtovsh)) then |
251 |
inb(i)=MIN(inb(i),k) |
252 |
endif |
253 |
535 continue |
254 |
530 continue |
255 |
|
256 |
!-- end convect3 |
257 |
|
258 |
! ori do 510 i=1,ncum |
259 |
! ori cape(i)=0.0 |
260 |
! ori capem(i)=0.0 |
261 |
! ori inb(i)=icb(i)+1 |
262 |
! ori inb1(i)=inb(i) |
263 |
! ori 510 continue |
264 |
! |
265 |
! Originial Code |
266 |
! |
267 |
! do 530 k=minorig+1,nl-1 |
268 |
! do 520 i=1,ncum |
269 |
! if(k.ge.(icb(i)+1))then |
270 |
! by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k) |
271 |
! byp=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1) |
272 |
! cape(i)=cape(i)+by |
273 |
! if(by.ge.0.0)inb1(i)=k+1 |
274 |
! if(cape(i).gt.0.0)then |
275 |
! inb(i)=k+1 |
276 |
! capem(i)=cape(i) |
277 |
! endif |
278 |
! endif |
279 |
!520 continue |
280 |
!530 continue |
281 |
! do 540 i=1,ncum |
282 |
! byp=(tvp(i,nl)-tv(i,nl))*dph(i,nl)/p(i,nl) |
283 |
! cape(i)=capem(i)+byp |
284 |
! defrac=capem(i)-cape(i) |
285 |
! defrac=max(defrac,0.001) |
286 |
! frac(i)=-cape(i)/defrac |
287 |
! frac(i)=min(frac(i),1.0) |
288 |
! frac(i)=max(frac(i),0.0) |
289 |
!540 continue |
290 |
! |
291 |
! K Emanuel fix |
292 |
! |
293 |
! call zilch(byp,ncum) |
294 |
! do 530 k=minorig+1,nl-1 |
295 |
! do 520 i=1,ncum |
296 |
! if(k.ge.(icb(i)+1))then |
297 |
! by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k) |
298 |
! cape(i)=cape(i)+by |
299 |
! if(by.ge.0.0)inb1(i)=k+1 |
300 |
! if(cape(i).gt.0.0)then |
301 |
! inb(i)=k+1 |
302 |
! capem(i)=cape(i) |
303 |
! byp(i)=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1) |
304 |
! endif |
305 |
! endif |
306 |
!520 continue |
307 |
!530 continue |
308 |
! do 540 i=1,ncum |
309 |
! inb(i)=max(inb(i),inb1(i)) |
310 |
! cape(i)=capem(i)+byp(i) |
311 |
! defrac=capem(i)-cape(i) |
312 |
! defrac=max(defrac,0.001) |
313 |
! frac(i)=-cape(i)/defrac |
314 |
! frac(i)=min(frac(i),1.0) |
315 |
! frac(i)=max(frac(i),0.0) |
316 |
!540 continue |
317 |
! |
318 |
! J Teixeira fix |
319 |
! |
320 |
! ori call zilch(byp,ncum) |
321 |
! ori do 515 i=1,ncum |
322 |
! ori lcape(i)=.true. |
323 |
! ori 515 continue |
324 |
! ori do 530 k=minorig+1,nl-1 |
325 |
! ori do 520 i=1,ncum |
326 |
! ori if(cape(i).lt.0.0)lcape(i)=.false. |
327 |
! ori if((k.ge.(icb(i)+1)).and.lcape(i))then |
328 |
! ori by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k) |
329 |
! ori byp(i)=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1) |
330 |
! ori cape(i)=cape(i)+by |
331 |
! ori if(by.ge.0.0)inb1(i)=k+1 |
332 |
! ori if(cape(i).gt.0.0)then |
333 |
! ori inb(i)=k+1 |
334 |
! ori capem(i)=cape(i) |
335 |
! ori endif |
336 |
! ori endif |
337 |
! ori 520 continue |
338 |
! ori 530 continue |
339 |
! ori do 540 i=1,ncum |
340 |
! ori cape(i)=capem(i)+byp(i) |
341 |
! ori defrac=capem(i)-cape(i) |
342 |
! ori defrac=max(defrac,0.001) |
343 |
! ori frac(i)=-cape(i)/defrac |
344 |
! ori frac(i)=min(frac(i),1.0) |
345 |
! ori frac(i)=max(frac(i),0.0) |
346 |
! ori 540 continue |
347 |
! |
348 |
!===================================================================== |
349 |
! --- CALCULATE LIQUID WATER STATIC ENERGY OF LIFTED PARCEL |
350 |
!===================================================================== |
351 |
! |
352 |
!ym do i=1,ncum*nlp |
353 |
!ym hp(i,1)=h(i,1) |
354 |
!ym enddo |
355 |
|
356 |
do k=1,nlp |
357 |
do i=1,ncum |
358 |
hp(i,k)=h(i,k) |
359 |
enddo |
360 |
enddo |
361 |
|
362 |
do 600 k=minorig+1,nl |
363 |
do 590 i=1,ncum |
364 |
if((k.ge.icb(i)).and.(k.le.inb(i)))then |
365 |
hp(i,k)=h(i,nk(i))+(lv(i,k)+(cpd-cpv)*t(i,k))*ep(i,k)*clw(i,k) |
366 |
endif |
367 |
590 continue |
368 |
600 continue |
369 |
|
370 |
return |
371 |
end |