1 |
guez |
3 |
! |
2 |
|
|
! $Header: /home/cvsroot/LMDZ4/libf/phylmd/cv_driver.F,v 1.3 2005/04/15 12:36:17 lmdzadmin Exp $ |
3 |
|
|
! |
4 |
|
|
SUBROUTINE cv_driver(len,nd,ndp1,ntra,iflag_con, |
5 |
|
|
& t1,q1,qs1,u1,v1,tra1, |
6 |
|
|
& p1,ph1,iflag1,ft1,fq1,fu1,fv1,ftra1, |
7 |
|
|
& precip1,VPrecip1, |
8 |
|
|
& cbmf1,sig1,w01, |
9 |
|
|
& icb1,inb1, |
10 |
|
|
& delt,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1, |
11 |
|
|
& da1,phi1,mp1) |
12 |
|
|
C |
13 |
|
|
use dimens_m |
14 |
|
|
use dimphy |
15 |
|
|
implicit none |
16 |
|
|
C |
17 |
|
|
C.............................START PROLOGUE............................ |
18 |
|
|
C |
19 |
|
|
C PARAMETERS: |
20 |
|
|
C Name Type Usage Description |
21 |
|
|
C ---------- ---------- ------- ---------------------------- |
22 |
|
|
C |
23 |
|
|
C len Integer Input first (i) dimension |
24 |
|
|
C nd Integer Input vertical (k) dimension |
25 |
|
|
C ndp1 Integer Input nd + 1 |
26 |
|
|
C ntra Integer Input number of tracors |
27 |
|
|
C iflag_con Integer Input version of convect (3/4) |
28 |
|
|
C t1 Real Input temperature |
29 |
|
|
C q1 Real Input specific hum |
30 |
|
|
C qs1 Real Input sat specific hum |
31 |
|
|
C u1 Real Input u-wind |
32 |
|
|
C v1 Real Input v-wind |
33 |
|
|
C tra1 Real Input tracors |
34 |
|
|
C p1 Real Input full level pressure |
35 |
|
|
C ph1 Real Input half level pressure |
36 |
|
|
C iflag1 Integer Output flag for Emanuel conditions |
37 |
|
|
C ft1 Real Output temp tend |
38 |
|
|
C fq1 Real Output spec hum tend |
39 |
|
|
C fu1 Real Output u-wind tend |
40 |
|
|
C fv1 Real Output v-wind tend |
41 |
|
|
C ftra1 Real Output tracor tend |
42 |
|
|
C precip1 Real Output precipitation |
43 |
|
|
C VPrecip1 Real Output vertical profile of precipitations |
44 |
|
|
C cbmf1 Real Output cloud base mass flux |
45 |
|
|
C sig1 Real In/Out section adiabatic updraft |
46 |
|
|
C w01 Real In/Out vertical velocity within adiab updraft |
47 |
|
|
C delt Real Input time step |
48 |
|
|
C Ma1 Real Output mass flux adiabatic updraft |
49 |
|
|
C upwd1 Real Output total upward mass flux (adiab+mixed) |
50 |
|
|
C dnwd1 Real Output saturated downward mass flux (mixed) |
51 |
|
|
C dnwd01 Real Output unsaturated downward mass flux |
52 |
|
|
C qcondc1 Real Output in-cld mixing ratio of condensed water |
53 |
|
|
C wd1 Real Output downdraft velocity scale for sfc fluxes |
54 |
|
|
C cape1 Real Output CAPE |
55 |
|
|
C |
56 |
|
|
C S. Bony, Mar 2002: |
57 |
|
|
C * Several modules corresponding to different physical processes |
58 |
|
|
C * Several versions of convect may be used: |
59 |
|
|
C - iflag_con=3: version lmd (previously named convect3) |
60 |
|
|
C - iflag_con=4: version 4.3b (vect. version, previously convect1/2) |
61 |
|
|
C + tard: - iflag_con=5: version lmd with ice (previously named convectg) |
62 |
|
|
C S. Bony, Oct 2002: |
63 |
|
|
C * Vectorization of convect3 (ie version lmd) |
64 |
|
|
C |
65 |
|
|
C..............................END PROLOGUE............................. |
66 |
|
|
c |
67 |
|
|
c |
68 |
|
|
|
69 |
|
|
integer len |
70 |
|
|
integer nd |
71 |
|
|
integer ndp1 |
72 |
|
|
integer noff |
73 |
guez |
12 |
integer, intent(in):: iflag_con |
74 |
guez |
3 |
integer ntra |
75 |
|
|
real t1(len,nd) |
76 |
|
|
real q1(len,nd) |
77 |
|
|
real qs1(len,nd) |
78 |
|
|
real u1(len,nd) |
79 |
|
|
real v1(len,nd) |
80 |
|
|
real p1(len,nd) |
81 |
|
|
real ph1(len,ndp1) |
82 |
|
|
integer iflag1(len) |
83 |
|
|
real ft1(len,nd) |
84 |
|
|
real fq1(len,nd) |
85 |
|
|
real fu1(len,nd) |
86 |
|
|
real fv1(len,nd) |
87 |
|
|
real precip1(len) |
88 |
|
|
real cbmf1(len) |
89 |
|
|
real VPrecip1(len,nd+1) |
90 |
|
|
real Ma1(len,nd) |
91 |
|
|
real upwd1(len,nd) |
92 |
|
|
real dnwd1(len,nd) |
93 |
|
|
real dnwd01(len,nd) |
94 |
|
|
|
95 |
|
|
real qcondc1(len,nd) ! cld |
96 |
|
|
real wd1(len) ! gust |
97 |
|
|
real cape1(len) |
98 |
|
|
|
99 |
|
|
real da1(len,nd),phi1(len,nd,nd),mp1(len,nd) |
100 |
|
|
real da(len,nd),phi(len,nd,nd),mp(len,nd) |
101 |
|
|
real, intent(in):: tra1(len,nd,ntra) |
102 |
|
|
real ftra1(len,nd,ntra) |
103 |
|
|
|
104 |
guez |
12 |
real, intent(in):: delt |
105 |
guez |
3 |
|
106 |
|
|
!------------------------------------------------------------------- |
107 |
|
|
! --- ARGUMENTS |
108 |
|
|
!------------------------------------------------------------------- |
109 |
|
|
! --- On input: |
110 |
|
|
! |
111 |
|
|
! t: Array of absolute temperature (K) of dimension ND, with first |
112 |
|
|
! index corresponding to lowest model level. Note that this array |
113 |
|
|
! will be altered by the subroutine if dry convective adjustment |
114 |
|
|
! occurs and if IPBL is not equal to 0. |
115 |
|
|
! |
116 |
|
|
! q: Array of specific humidity (gm/gm) of dimension ND, with first |
117 |
|
|
! index corresponding to lowest model level. Must be defined |
118 |
|
|
! at same grid levels as T. Note that this array will be altered |
119 |
|
|
! if dry convective adjustment occurs and if IPBL is not equal to 0. |
120 |
|
|
! |
121 |
|
|
! qs: Array of saturation specific humidity of dimension ND, with first |
122 |
|
|
! index corresponding to lowest model level. Must be defined |
123 |
|
|
! at same grid levels as T. Note that this array will be altered |
124 |
|
|
! if dry convective adjustment occurs and if IPBL is not equal to 0. |
125 |
|
|
! |
126 |
|
|
! u: Array of zonal wind velocity (m/s) of dimension ND, witth first |
127 |
|
|
! index corresponding with the lowest model level. Defined at |
128 |
|
|
! same levels as T. Note that this array will be altered if |
129 |
|
|
! dry convective adjustment occurs and if IPBL is not equal to 0. |
130 |
|
|
! |
131 |
|
|
! v: Same as u but for meridional velocity. |
132 |
|
|
! |
133 |
|
|
! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA), |
134 |
|
|
! where NTRA is the number of different tracers. If no |
135 |
|
|
! convective tracer transport is needed, define a dummy |
136 |
|
|
! input array of dimension (ND,1). Tracers are defined at |
137 |
|
|
! same vertical levels as T. Note that this array will be altered |
138 |
|
|
! if dry convective adjustment occurs and if IPBL is not equal to 0. |
139 |
|
|
! |
140 |
|
|
! p: Array of pressure (mb) of dimension ND, with first |
141 |
|
|
! index corresponding to lowest model level. Must be defined |
142 |
|
|
! at same grid levels as T. |
143 |
|
|
! |
144 |
|
|
! ph: Array of pressure (mb) of dimension ND+1, with first index |
145 |
|
|
! corresponding to lowest level. These pressures are defined at |
146 |
|
|
! levels intermediate between those of P, T, Q and QS. The first |
147 |
|
|
! value of PH should be greater than (i.e. at a lower level than) |
148 |
|
|
! the first value of the array P. |
149 |
|
|
! |
150 |
|
|
! nl: The maximum number of levels to which convection can penetrate, plus 1. |
151 |
|
|
! NL MUST be less than or equal to ND-1. |
152 |
|
|
! |
153 |
|
|
! delt: The model time step (sec) between calls to CONVECT |
154 |
|
|
! |
155 |
|
|
!---------------------------------------------------------------------------- |
156 |
|
|
! --- On Output: |
157 |
|
|
! |
158 |
|
|
! iflag: An output integer whose value denotes the following: |
159 |
|
|
! VALUE INTERPRETATION |
160 |
|
|
! ----- -------------- |
161 |
|
|
! 0 Moist convection occurs. |
162 |
|
|
! 1 Moist convection occurs, but a CFL condition |
163 |
|
|
! on the subsidence warming is violated. This |
164 |
|
|
! does not cause the scheme to terminate. |
165 |
|
|
! 2 Moist convection, but no precip because ep(inb) lt 0.0001 |
166 |
|
|
! 3 No moist convection because new cbmf is 0 and old cbmf is 0. |
167 |
|
|
! 4 No moist convection; atmosphere is not |
168 |
|
|
! unstable |
169 |
|
|
! 6 No moist convection because ihmin le minorig. |
170 |
|
|
! 7 No moist convection because unreasonable |
171 |
|
|
! parcel level temperature or specific humidity. |
172 |
|
|
! 8 No moist convection: lifted condensation |
173 |
|
|
! level is above the 200 mb level. |
174 |
|
|
! 9 No moist convection: cloud base is higher |
175 |
|
|
! then the level NL-1. |
176 |
|
|
! |
177 |
|
|
! ft: Array of temperature tendency (K/s) of dimension ND, defined at same |
178 |
|
|
! grid levels as T, Q, QS and P. |
179 |
|
|
! |
180 |
|
|
! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND, |
181 |
|
|
! defined at same grid levels as T, Q, QS and P. |
182 |
|
|
! |
183 |
|
|
! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND, |
184 |
|
|
! defined at same grid levels as T. |
185 |
|
|
! |
186 |
|
|
! fv: Same as FU, but for forcing of meridional velocity. |
187 |
|
|
! |
188 |
|
|
! ftra: Array of forcing of tracer content, in tracer mixing ratio per |
189 |
|
|
! second, defined at same levels as T. Dimensioned (ND,NTRA). |
190 |
|
|
! |
191 |
|
|
! precip: Scalar convective precipitation rate (mm/day). |
192 |
|
|
! |
193 |
|
|
! VPrecip: Vertical profile of convective precipitation (kg/m2/s). |
194 |
|
|
! |
195 |
|
|
! wd: A convective downdraft velocity scale. For use in surface |
196 |
|
|
! flux parameterizations. See convect.ps file for details. |
197 |
|
|
! |
198 |
|
|
! tprime: A convective downdraft temperature perturbation scale (K). |
199 |
|
|
! For use in surface flux parameterizations. See convect.ps |
200 |
|
|
! file for details. |
201 |
|
|
! |
202 |
|
|
! qprime: A convective downdraft specific humidity |
203 |
|
|
! perturbation scale (gm/gm). |
204 |
|
|
! For use in surface flux parameterizations. See convect.ps |
205 |
|
|
! file for details. |
206 |
|
|
! |
207 |
|
|
! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST |
208 |
|
|
! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT |
209 |
|
|
! ITS NEXT CALL. That is, the value of CBMF must be "remembered" |
210 |
|
|
! by the calling program between calls to CONVECT. |
211 |
|
|
! |
212 |
|
|
! det: Array of detrainment mass flux of dimension ND. |
213 |
|
|
! |
214 |
|
|
!------------------------------------------------------------------- |
215 |
|
|
c |
216 |
|
|
c Local arrays |
217 |
|
|
c |
218 |
|
|
|
219 |
|
|
integer i,k,n,il,j |
220 |
|
|
integer icbmax |
221 |
|
|
integer nk1(klon) |
222 |
|
|
integer icb1(klon) |
223 |
|
|
integer inb1(klon) |
224 |
|
|
integer icbs1(klon) |
225 |
|
|
|
226 |
|
|
real plcl1(klon) |
227 |
|
|
real tnk1(klon) |
228 |
|
|
real qnk1(klon) |
229 |
|
|
real gznk1(klon) |
230 |
|
|
real pnk1(klon) |
231 |
|
|
real qsnk1(klon) |
232 |
|
|
real pbase1(klon) |
233 |
|
|
real buoybase1(klon) |
234 |
|
|
|
235 |
|
|
real lv1(klon,klev) |
236 |
|
|
real cpn1(klon,klev) |
237 |
|
|
real tv1(klon,klev) |
238 |
|
|
real gz1(klon,klev) |
239 |
|
|
real hm1(klon,klev) |
240 |
|
|
real h1(klon,klev) |
241 |
|
|
real tp1(klon,klev) |
242 |
|
|
real tvp1(klon,klev) |
243 |
|
|
real clw1(klon,klev) |
244 |
|
|
real sig1(klon,klev) |
245 |
|
|
real w01(klon,klev) |
246 |
|
|
real th1(klon,klev) |
247 |
|
|
c |
248 |
|
|
integer ncum |
249 |
|
|
c |
250 |
|
|
c (local) compressed fields: |
251 |
|
|
c |
252 |
|
|
integer nloc |
253 |
|
|
parameter (nloc=klon) ! pour l'instant |
254 |
|
|
|
255 |
|
|
integer idcum(nloc) |
256 |
|
|
integer iflag(nloc),nk(nloc),icb(nloc) |
257 |
|
|
integer nent(nloc,klev) |
258 |
|
|
integer icbs(nloc) |
259 |
|
|
integer inb(nloc), inbis(nloc) |
260 |
|
|
|
261 |
|
|
real cbmf(nloc),plcl(nloc),tnk(nloc),qnk(nloc),gznk(nloc) |
262 |
|
|
real t(nloc,klev),q(nloc,klev),qs(nloc,klev) |
263 |
|
|
real u(nloc,klev),v(nloc,klev) |
264 |
|
|
real gz(nloc,klev),h(nloc,klev),lv(nloc,klev),cpn(nloc,klev) |
265 |
|
|
real p(nloc,klev),ph(nloc,klev+1),tv(nloc,klev),tp(nloc,klev) |
266 |
|
|
real clw(nloc,klev) |
267 |
|
|
real dph(nloc,klev) |
268 |
|
|
real pbase(nloc), buoybase(nloc), th(nloc,klev) |
269 |
|
|
real tvp(nloc,klev) |
270 |
|
|
real sig(nloc,klev), w0(nloc,klev) |
271 |
|
|
real hp(nloc,klev), ep(nloc,klev), sigp(nloc,klev) |
272 |
|
|
real frac(nloc), buoy(nloc,klev) |
273 |
|
|
real cape(nloc) |
274 |
|
|
real m(nloc,klev), ment(nloc,klev,klev), qent(nloc,klev,klev) |
275 |
|
|
real uent(nloc,klev,klev), vent(nloc,klev,klev) |
276 |
|
|
real ments(nloc,klev,klev), qents(nloc,klev,klev) |
277 |
|
|
real sij(nloc,klev,klev), elij(nloc,klev,klev) |
278 |
|
|
real qp(nloc,klev), up(nloc,klev), vp(nloc,klev) |
279 |
|
|
real wt(nloc,klev), water(nloc,klev), evap(nloc,klev) |
280 |
|
|
real b(nloc,klev), ft(nloc,klev), fq(nloc,klev) |
281 |
|
|
real fu(nloc,klev), fv(nloc,klev) |
282 |
|
|
real upwd(nloc,klev), dnwd(nloc,klev), dnwd0(nloc,klev) |
283 |
|
|
real Ma(nloc,klev), mike(nloc,klev), tls(nloc,klev) |
284 |
|
|
real tps(nloc,klev), qprime(nloc), tprime(nloc) |
285 |
|
|
real precip(nloc) |
286 |
|
|
real VPrecip(nloc,klev+1) |
287 |
|
|
real tra(nloc,klev,ntra), trap(nloc,klev,ntra) |
288 |
|
|
real ftra(nloc,klev,ntra), traent(nloc,klev,klev,ntra) |
289 |
|
|
real qcondc(nloc,klev) ! cld |
290 |
|
|
real wd(nloc) ! gust |
291 |
|
|
|
292 |
|
|
!------------------------------------------------------------------- |
293 |
|
|
! --- SET CONSTANTS AND PARAMETERS |
294 |
|
|
!------------------------------------------------------------------- |
295 |
|
|
|
296 |
|
|
c -- set simulation flags: |
297 |
|
|
c (common cvflag) |
298 |
|
|
|
299 |
|
|
CALL cv_flag |
300 |
|
|
|
301 |
|
|
c -- set thermodynamical constants: |
302 |
|
|
c (common cvthermo) |
303 |
|
|
|
304 |
|
|
CALL cv_thermo(iflag_con) |
305 |
|
|
|
306 |
|
|
c -- set convect parameters |
307 |
|
|
c |
308 |
|
|
c includes microphysical parameters and parameters that |
309 |
|
|
c control the rate of approach to quasi-equilibrium) |
310 |
|
|
c (common cvparam) |
311 |
|
|
|
312 |
|
|
if (iflag_con.eq.3) then |
313 |
|
|
CALL cv3_param(nd,delt) |
314 |
|
|
endif |
315 |
|
|
|
316 |
|
|
if (iflag_con.eq.4) then |
317 |
|
|
CALL cv_param(nd) |
318 |
|
|
endif |
319 |
|
|
|
320 |
|
|
!--------------------------------------------------------------------- |
321 |
|
|
! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS |
322 |
|
|
!--------------------------------------------------------------------- |
323 |
|
|
|
324 |
|
|
do 20 k=1,nd |
325 |
|
|
do 10 i=1,len |
326 |
|
|
ft1(i,k)=0.0 |
327 |
|
|
fq1(i,k)=0.0 |
328 |
|
|
fu1(i,k)=0.0 |
329 |
|
|
fv1(i,k)=0.0 |
330 |
|
|
tvp1(i,k)=0.0 |
331 |
|
|
tp1(i,k)=0.0 |
332 |
|
|
clw1(i,k)=0.0 |
333 |
|
|
cym |
334 |
|
|
clw(i,k)=0.0 |
335 |
|
|
gz1(i,k) = 0. |
336 |
|
|
VPrecip1(i,k) = 0. |
337 |
|
|
Ma1(i,k)=0.0 |
338 |
|
|
upwd1(i,k)=0.0 |
339 |
|
|
dnwd1(i,k)=0.0 |
340 |
|
|
dnwd01(i,k)=0.0 |
341 |
|
|
qcondc1(i,k)=0.0 |
342 |
|
|
10 continue |
343 |
|
|
20 continue |
344 |
|
|
|
345 |
|
|
do 30 j=1,ntra |
346 |
|
|
do 31 k=1,nd |
347 |
|
|
do 32 i=1,len |
348 |
|
|
ftra1(i,k,j)=0.0 |
349 |
|
|
32 continue |
350 |
|
|
31 continue |
351 |
|
|
30 continue |
352 |
|
|
|
353 |
|
|
do 60 i=1,len |
354 |
|
|
precip1(i)=0.0 |
355 |
|
|
iflag1(i)=0 |
356 |
|
|
wd1(i)=0.0 |
357 |
|
|
cape1(i)=0.0 |
358 |
|
|
VPrecip1(i,nd+1)=0.0 |
359 |
|
|
60 continue |
360 |
|
|
|
361 |
|
|
if (iflag_con.eq.3) then |
362 |
|
|
do il=1,len |
363 |
|
|
sig1(il,nd)=sig1(il,nd)+1. |
364 |
|
|
sig1(il,nd)=amin1(sig1(il,nd),12.1) |
365 |
|
|
enddo |
366 |
|
|
endif |
367 |
|
|
|
368 |
|
|
!-------------------------------------------------------------------- |
369 |
|
|
! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY |
370 |
|
|
!-------------------------------------------------------------------- |
371 |
|
|
|
372 |
|
|
if (iflag_con.eq.3) then |
373 |
|
|
CALL cv3_prelim(len,nd,ndp1,t1,q1,p1,ph1 ! nd->na |
374 |
|
|
o ,lv1,cpn1,tv1,gz1,h1,hm1,th1) |
375 |
|
|
endif |
376 |
|
|
|
377 |
|
|
if (iflag_con.eq.4) then |
378 |
|
|
CALL cv_prelim(len,nd,ndp1,t1,q1,p1,ph1 |
379 |
|
|
o ,lv1,cpn1,tv1,gz1,h1,hm1) |
380 |
|
|
endif |
381 |
|
|
|
382 |
|
|
!-------------------------------------------------------------------- |
383 |
|
|
! --- CONVECTIVE FEED |
384 |
|
|
!-------------------------------------------------------------------- |
385 |
|
|
|
386 |
|
|
if (iflag_con.eq.3) then |
387 |
|
|
CALL cv3_feed(len,nd,t1,q1,qs1,p1,ph1,hm1,gz1 ! nd->na |
388 |
|
|
o ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1) |
389 |
|
|
endif |
390 |
|
|
|
391 |
|
|
if (iflag_con.eq.4) then |
392 |
|
|
CALL cv_feed(len,nd,t1,q1,qs1,p1,hm1,gz1 |
393 |
|
|
o ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1) |
394 |
|
|
endif |
395 |
|
|
|
396 |
|
|
!-------------------------------------------------------------------- |
397 |
|
|
! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part |
398 |
|
|
! (up through ICB for convect4, up through ICB+1 for convect3) |
399 |
|
|
! Calculates the lifted parcel virtual temperature at nk, the |
400 |
|
|
! actual temperature, and the adiabatic liquid water content. |
401 |
|
|
!-------------------------------------------------------------------- |
402 |
|
|
|
403 |
|
|
if (iflag_con.eq.3) then |
404 |
|
|
CALL cv3_undilute1(len,nd,t1,q1,qs1,gz1,plcl1,p1,nk1,icb1 ! nd->na |
405 |
|
|
o ,tp1,tvp1,clw1,icbs1) |
406 |
|
|
endif |
407 |
|
|
|
408 |
|
|
if (iflag_con.eq.4) then |
409 |
|
|
CALL cv_undilute1(len,nd,t1,q1,qs1,gz1,p1,nk1,icb1,icbmax |
410 |
|
|
: ,tp1,tvp1,clw1) |
411 |
|
|
endif |
412 |
|
|
|
413 |
|
|
!------------------------------------------------------------------- |
414 |
|
|
! --- TRIGGERING |
415 |
|
|
!------------------------------------------------------------------- |
416 |
|
|
|
417 |
|
|
if (iflag_con.eq.3) then |
418 |
|
|
CALL cv3_trigger(len,nd,icb1,plcl1,p1,th1,tv1,tvp1 ! nd->na |
419 |
|
|
o ,pbase1,buoybase1,iflag1,sig1,w01) |
420 |
|
|
endif |
421 |
|
|
|
422 |
|
|
if (iflag_con.eq.4) then |
423 |
|
|
CALL cv_trigger(len,nd,icb1,cbmf1,tv1,tvp1,iflag1) |
424 |
|
|
endif |
425 |
|
|
|
426 |
|
|
!===================================================================== |
427 |
|
|
! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY |
428 |
|
|
!===================================================================== |
429 |
|
|
|
430 |
|
|
ncum=0 |
431 |
|
|
do 400 i=1,len |
432 |
|
|
if(iflag1(i).eq.0)then |
433 |
|
|
ncum=ncum+1 |
434 |
|
|
idcum(ncum)=i |
435 |
|
|
endif |
436 |
|
|
400 continue |
437 |
|
|
|
438 |
|
|
c print*,'klon, ncum = ',len,ncum |
439 |
|
|
|
440 |
|
|
IF (ncum.gt.0) THEN |
441 |
|
|
|
442 |
|
|
!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
443 |
|
|
! --- COMPRESS THE FIELDS |
444 |
|
|
! (-> vectorization over convective gridpoints) |
445 |
|
|
!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
446 |
|
|
|
447 |
|
|
if (iflag_con.eq.3) then |
448 |
|
|
CALL cv3_compress( len,nloc,ncum,nd,ntra |
449 |
|
|
: ,iflag1,nk1,icb1,icbs1 |
450 |
|
|
: ,plcl1,tnk1,qnk1,gznk1,pbase1,buoybase1 |
451 |
|
|
: ,t1,q1,qs1,u1,v1,gz1,th1 |
452 |
|
|
: ,tra1 |
453 |
|
|
: ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1 |
454 |
|
|
: ,sig1,w01 |
455 |
|
|
o ,iflag,nk,icb,icbs |
456 |
|
|
o ,plcl,tnk,qnk,gznk,pbase,buoybase |
457 |
|
|
o ,t,q,qs,u,v,gz,th |
458 |
|
|
o ,tra |
459 |
|
|
o ,h,lv,cpn,p,ph,tv,tp,tvp,clw |
460 |
|
|
o ,sig,w0 ) |
461 |
|
|
endif |
462 |
|
|
|
463 |
|
|
if (iflag_con.eq.4) then |
464 |
|
|
CALL cv_compress( len,nloc,ncum,nd |
465 |
|
|
: ,iflag1,nk1,icb1 |
466 |
|
|
: ,cbmf1,plcl1,tnk1,qnk1,gznk1 |
467 |
|
|
: ,t1,q1,qs1,u1,v1,gz1 |
468 |
|
|
: ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1 |
469 |
|
|
o ,iflag,nk,icb |
470 |
|
|
o ,cbmf,plcl,tnk,qnk,gznk |
471 |
|
|
o ,t,q,qs,u,v,gz,h,lv,cpn,p,ph,tv,tp,tvp,clw |
472 |
|
|
o ,dph ) |
473 |
|
|
endif |
474 |
|
|
|
475 |
|
|
!------------------------------------------------------------------- |
476 |
|
|
! --- UNDILUTE (ADIABATIC) UPDRAFT / second part : |
477 |
|
|
! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
478 |
|
|
! --- & |
479 |
|
|
! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE |
480 |
|
|
! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD |
481 |
|
|
! --- & |
482 |
|
|
! --- FIND THE LEVEL OF NEUTRAL BUOYANCY |
483 |
|
|
!------------------------------------------------------------------- |
484 |
|
|
|
485 |
|
|
if (iflag_con.eq.3) then |
486 |
|
|
CALL cv3_undilute2(nloc,ncum,nd,icb,icbs,nk !na->nd |
487 |
|
|
: ,tnk,qnk,gznk,t,q,qs,gz |
488 |
|
|
: ,p,h,tv,lv,pbase,buoybase,plcl |
489 |
|
|
o ,inb,tp,tvp,clw,hp,ep,sigp,buoy) |
490 |
|
|
endif |
491 |
|
|
|
492 |
|
|
if (iflag_con.eq.4) then |
493 |
|
|
CALL cv_undilute2(nloc,ncum,nd,icb,nk |
494 |
|
|
: ,tnk,qnk,gznk,t,q,qs,gz |
495 |
|
|
: ,p,dph,h,tv,lv |
496 |
|
|
o ,inb,inbis,tp,tvp,clw,hp,ep,sigp,frac) |
497 |
|
|
endif |
498 |
|
|
|
499 |
|
|
!------------------------------------------------------------------- |
500 |
|
|
! --- CLOSURE |
501 |
|
|
!------------------------------------------------------------------- |
502 |
|
|
|
503 |
|
|
if (iflag_con.eq.3) then |
504 |
|
|
CALL cv3_closure(nloc,ncum,nd,icb,inb ! na->nd |
505 |
|
|
: ,pbase,p,ph,tv,buoy |
506 |
|
|
o ,sig,w0,cape,m) |
507 |
|
|
endif |
508 |
|
|
|
509 |
|
|
if (iflag_con.eq.4) then |
510 |
|
|
CALL cv_closure(nloc,ncum,nd,nk,icb |
511 |
|
|
: ,tv,tvp,p,ph,dph,plcl,cpn |
512 |
|
|
o ,iflag,cbmf) |
513 |
|
|
endif |
514 |
|
|
|
515 |
|
|
!------------------------------------------------------------------- |
516 |
|
|
! --- MIXING |
517 |
|
|
!------------------------------------------------------------------- |
518 |
|
|
|
519 |
|
|
if (iflag_con.eq.3) then |
520 |
|
|
CALL cv3_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb ! na->nd |
521 |
|
|
: ,ph,t,q,qs,u,v,tra,h,lv,qnk |
522 |
|
|
: ,hp,tv,tvp,ep,clw,m,sig |
523 |
|
|
o ,ment,qent,uent,vent,sij,elij,ments,qents,traent) |
524 |
|
|
endif |
525 |
|
|
|
526 |
|
|
if (iflag_con.eq.4) then |
527 |
|
|
CALL cv_mixing(nloc,ncum,nd,icb,nk,inb,inbis |
528 |
|
|
: ,ph,t,q,qs,u,v,h,lv,qnk |
529 |
|
|
: ,hp,tv,tvp,ep,clw,cbmf |
530 |
|
|
o ,m,ment,qent,uent,vent,nent,sij,elij) |
531 |
|
|
endif |
532 |
|
|
|
533 |
|
|
!------------------------------------------------------------------- |
534 |
|
|
! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS |
535 |
|
|
!------------------------------------------------------------------- |
536 |
|
|
|
537 |
|
|
if (iflag_con.eq.3) then |
538 |
|
|
CALL cv3_unsat(nloc,ncum,nd,nd,ntra,icb,inb ! na->nd |
539 |
|
|
: ,t,q,qs,gz,u,v,tra,p,ph |
540 |
|
|
: ,th,tv,lv,cpn,ep,sigp,clw |
541 |
|
|
: ,m,ment,elij,delt,plcl |
542 |
|
|
o ,mp,qp,up,vp,trap,wt,water,evap,b) |
543 |
|
|
endif |
544 |
|
|
|
545 |
|
|
if (iflag_con.eq.4) then |
546 |
|
|
CALL cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph |
547 |
|
|
: ,h,lv,ep,sigp,clw,m,ment,elij |
548 |
|
|
o ,iflag,mp,qp,up,vp,wt,water,evap) |
549 |
|
|
endif |
550 |
|
|
|
551 |
|
|
!------------------------------------------------------------------- |
552 |
|
|
! --- YIELD |
553 |
|
|
! (tendencies, precipitation, variables of interface with other |
554 |
|
|
! processes, etc) |
555 |
|
|
!------------------------------------------------------------------- |
556 |
|
|
|
557 |
|
|
if (iflag_con.eq.3) then |
558 |
|
|
CALL cv3_yield(nloc,ncum,nd,nd,ntra ! na->nd |
559 |
|
|
: ,icb,inb,delt |
560 |
|
|
: ,t,q,u,v,tra,gz,p,ph,h,hp,lv,cpn,th |
561 |
|
|
: ,ep,clw,m,tp,mp,qp,up,vp,trap |
562 |
|
|
: ,wt,water,evap,b |
563 |
|
|
: ,ment,qent,uent,vent,nent,elij,traent,sig |
564 |
|
|
: ,tv,tvp |
565 |
|
|
o ,iflag,precip,VPrecip,ft,fq,fu,fv,ftra |
566 |
|
|
o ,upwd,dnwd,dnwd0,ma,mike,tls,tps,qcondc,wd) |
567 |
|
|
endif |
568 |
|
|
|
569 |
|
|
if (iflag_con.eq.4) then |
570 |
|
|
CALL cv_yield(nloc,ncum,nd,nk,icb,inb,delt |
571 |
|
|
: ,t,q,u,v,gz,p,ph,h,hp,lv,cpn |
572 |
|
|
: ,ep,clw,frac,m,mp,qp,up,vp |
573 |
|
|
: ,wt,water,evap |
574 |
|
|
: ,ment,qent,uent,vent,nent,elij |
575 |
|
|
: ,tv,tvp |
576 |
|
|
o ,iflag,wd,qprime,tprime |
577 |
|
|
o ,precip,cbmf,ft,fq,fu,fv,Ma,qcondc) |
578 |
|
|
endif |
579 |
|
|
|
580 |
|
|
!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
581 |
|
|
! --- passive tracers |
582 |
|
|
!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
583 |
|
|
|
584 |
|
|
if (iflag_con.eq.3) then |
585 |
|
|
CALL cv3_tracer(nloc,len,ncum,nd,nd, |
586 |
|
|
: ment,sij,da,phi) |
587 |
|
|
endif |
588 |
|
|
|
589 |
|
|
!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
590 |
|
|
! --- UNCOMPRESS THE FIELDS |
591 |
|
|
!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
592 |
|
|
c set iflag1 =42 for non convective points |
593 |
|
|
do i=1,len |
594 |
|
|
iflag1(i)=42 |
595 |
|
|
end do |
596 |
|
|
c |
597 |
|
|
if (iflag_con.eq.3) then |
598 |
|
|
CALL cv3_uncompress(nloc,len,ncum,nd,ntra,idcum |
599 |
|
|
: ,iflag |
600 |
|
|
: ,precip,VPrecip,sig,w0 |
601 |
|
|
: ,ft,fq,fu,fv,ftra |
602 |
|
|
: ,inb |
603 |
|
|
: ,Ma,upwd,dnwd,dnwd0,qcondc,wd,cape |
604 |
|
|
: ,da,phi,mp |
605 |
|
|
o ,iflag1 |
606 |
|
|
o ,precip1,VPrecip1,sig1,w01 |
607 |
|
|
o ,ft1,fq1,fu1,fv1,ftra1 |
608 |
|
|
o ,inb1 |
609 |
|
|
o ,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1 |
610 |
|
|
o ,da1,phi1,mp1) |
611 |
|
|
endif |
612 |
|
|
|
613 |
|
|
if (iflag_con.eq.4) then |
614 |
|
|
CALL cv_uncompress(nloc,len,ncum,nd,idcum |
615 |
|
|
: ,iflag |
616 |
|
|
: ,precip,cbmf |
617 |
|
|
: ,ft,fq,fu,fv |
618 |
|
|
: ,Ma,qcondc |
619 |
|
|
o ,iflag1 |
620 |
|
|
o ,precip1,cbmf1 |
621 |
|
|
o ,ft1,fq1,fu1,fv1 |
622 |
|
|
o ,Ma1,qcondc1 ) |
623 |
|
|
endif |
624 |
|
|
|
625 |
|
|
ENDIF ! ncum>0 |
626 |
|
|
|
627 |
|
|
9999 continue |
628 |
|
|
|
629 |
|
|
return |
630 |
|
|
end |
631 |
|
|
|
632 |
|
|
!================================================================== |
633 |
|
|
SUBROUTINE cv_flag |
634 |
|
|
implicit none |
635 |
|
|
|
636 |
|
|
include "cvflag.h" |
637 |
|
|
|
638 |
|
|
c -- si .TRUE., on rend la gravite plus explicite et eventuellement |
639 |
|
|
c differente de 10.0 dans convect3: |
640 |
|
|
cvflag_grav = .TRUE. |
641 |
|
|
|
642 |
|
|
return |
643 |
|
|
end |
644 |
|
|
|
645 |
|
|
!================================================================== |
646 |
|
|
SUBROUTINE cv_thermo(iflag_con) |
647 |
|
|
use YOMCST |
648 |
|
|
implicit none |
649 |
|
|
|
650 |
|
|
c------------------------------------------------------------- |
651 |
|
|
c Set thermodynamical constants for convectL |
652 |
|
|
c------------------------------------------------------------- |
653 |
|
|
|
654 |
|
|
include "cvthermo.h" |
655 |
|
|
|
656 |
guez |
12 |
integer, intent(in):: iflag_con |
657 |
guez |
3 |
|
658 |
|
|
|
659 |
|
|
c original set from convect: |
660 |
|
|
if (iflag_con.eq.4) then |
661 |
|
|
cpd=1005.7 |
662 |
|
|
cpv=1870.0 |
663 |
|
|
cl=4190.0 |
664 |
|
|
rrv=461.5 |
665 |
|
|
rrd=287.04 |
666 |
|
|
lv0=2.501E6 |
667 |
|
|
g=9.8 |
668 |
|
|
t0=273.15 |
669 |
|
|
grav=g |
670 |
|
|
endif |
671 |
|
|
|
672 |
|
|
c constants consistent with LMDZ: |
673 |
|
|
if (iflag_con.eq.3) then |
674 |
|
|
cpd = RCPD |
675 |
|
|
cpv = RCPV |
676 |
|
|
cl = RCW |
677 |
|
|
rrv = RV |
678 |
|
|
rrd = RD |
679 |
|
|
lv0 = RLVTT |
680 |
|
|
g = RG ! not used in convect3 |
681 |
|
|
c ori t0 = RTT |
682 |
|
|
t0 = 273.15 ! convect3 (RTT=273.16) |
683 |
|
|
c maf grav= 10. ! implicitely or explicitely used in convect3 |
684 |
|
|
grav= g ! implicitely or explicitely used in convect3 |
685 |
|
|
endif |
686 |
|
|
|
687 |
|
|
rowl=1000.0 !(a quelle variable de YOMCST cela correspond-il?) |
688 |
|
|
|
689 |
|
|
clmcpv=cl-cpv |
690 |
|
|
clmcpd=cl-cpd |
691 |
|
|
cpdmcp=cpd-cpv |
692 |
|
|
cpvmcpd=cpv-cpd |
693 |
|
|
cpvmcl=cl-cpv ! for convect3 |
694 |
|
|
eps=rrd/rrv |
695 |
|
|
epsi=1.0/eps |
696 |
|
|
epsim1=epsi-1.0 |
697 |
|
|
c ginv=1.0/g |
698 |
|
|
ginv=1.0/grav |
699 |
|
|
hrd=0.5*rrd |
700 |
|
|
|
701 |
|
|
return |
702 |
|
|
end |
703 |
|
|
|