1 | |
---|
2 | ! $Id: cva_driver.F90 2925 2017-06-30 10:17:42Z lguez $ |
---|
3 | |
---|
4 | SUBROUTINE cva_driver(len, nd, ndp1, ntra, nloc, k_upper, & |
---|
5 | iflag_con, iflag_mix, iflag_ice_thermo, iflag_clos, ok_conserv_q, & |
---|
6 | !! delt, t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, & ! jyg |
---|
7 | delt, comp_threshold, & ! jyg |
---|
8 | t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, & ! jyg |
---|
9 | u1, v1, tra1, & |
---|
10 | p1, ph1, & |
---|
11 | Ale1, Alp1, omega1, & |
---|
12 | sig1feed1, sig2feed1, wght1, & |
---|
13 | iflag1, ft1, fq1, fu1, fv1, ftra1, & |
---|
14 | precip1, kbas1, ktop1, & |
---|
15 | cbmf1, plcl1, plfc1, wbeff1, & |
---|
16 | sig1, w01, & !input/output |
---|
17 | ptop21, sigd1, & |
---|
18 | ma1, mip1, Vprecip1, Vprecipi1, upwd1, dnwd1, dnwd01, & ! jyg |
---|
19 | qcondc1, wd1, & |
---|
20 | cape1, cin1, tvp1, & |
---|
21 | ftd1, fqd1, & |
---|
22 | Plim11, Plim21, asupmax1, supmax01, asupmaxmin1, & |
---|
23 | lalim_conv1, & |
---|
24 | !! da1,phi1,mp1,phi21,d1a1,dam1,sigij1,clw1, & ! RomP |
---|
25 | !! elij1,evap1,ep1,epmlmMm1,eplaMm1, & ! RomP |
---|
26 | da1, phi1, mp1, phi21, d1a1, dam1, sigij1, wghti1, & ! RomP, RL |
---|
27 | clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, & ! RomP, RL |
---|
28 | wdtrainA1, wdtrainM1, qtc1, sigt1, tau_cld_cv, & |
---|
29 | coefw_cld_cv, & ! RomP, AJ |
---|
30 | epmax_diag1) ! epmax_cape |
---|
31 | ! ************************************************************** |
---|
32 | ! * |
---|
33 | ! CV_DRIVER * |
---|
34 | ! * |
---|
35 | ! * |
---|
36 | ! written by : Sandrine Bony-Lena , 17/05/2003, 11.19.41 * |
---|
37 | ! modified by : * |
---|
38 | ! ************************************************************** |
---|
39 | ! ************************************************************** |
---|
40 | |
---|
41 | USE print_control_mod, ONLY: prt_level, lunout |
---|
42 | USE add_phys_tend_mod, ONLY: fl_cor_ebil |
---|
43 | IMPLICIT NONE |
---|
44 | |
---|
45 | ! .............................START PROLOGUE............................ |
---|
46 | |
---|
47 | |
---|
48 | ! All argument names (except len,nd,ntra,nloc,delt and the flags) have a "1" appended. |
---|
49 | ! The "1" is removed for the corresponding compressed variables. |
---|
50 | ! PARAMETERS: |
---|
51 | ! Name Type Usage Description |
---|
52 | ! ---------- ---------- ------- ---------------------------- |
---|
53 | |
---|
54 | ! len Integer Input first (i) dimension |
---|
55 | ! nd Integer Input vertical (k) dimension |
---|
56 | ! ndp1 Integer Input nd + 1 |
---|
57 | ! ntra Integer Input number of tracors |
---|
58 | ! nloc Integer Input dimension of arrays for compressed fields |
---|
59 | ! k_upper Integer Input upmost level for vertical loops |
---|
60 | ! iflag_con Integer Input version of convect (3/4) |
---|
61 | ! iflag_mix Integer Input version of mixing (0/1/2) |
---|
62 | ! iflag_ice_thermo Integer Input accounting for ice thermodynamics (0/1) |
---|
63 | ! iflag_clos Integer Input version of closure (0/1) |
---|
64 | ! tau_cld_cv Real Input characteristic time of dissipation of mixing fluxes |
---|
65 | ! coefw_cld_cv Real Input coefficient for updraft velocity in convection |
---|
66 | ! ok_conserv_q Logical Input when true corrections for water conservation are swtiched on |
---|
67 | ! delt Real Input time step |
---|
68 | ! comp_threshold Real Input threshold on the fraction of convective points below which |
---|
69 | ! fields are compressed |
---|
70 | ! t1 Real Input temperature (sat draught envt) |
---|
71 | ! q1 Real Input specific hum (sat draught envt) |
---|
72 | ! qs1 Real Input sat specific hum (sat draught envt) |
---|
73 | ! t1_wake Real Input temperature (unsat draught envt) |
---|
74 | ! q1_wake Real Input specific hum(unsat draught envt) |
---|
75 | ! qs1_wake Real Input sat specific hum(unsat draughts envt) |
---|
76 | ! s1_wake Real Input fractionnal area covered by wakes |
---|
77 | ! u1 Real Input u-wind |
---|
78 | ! v1 Real Input v-wind |
---|
79 | ! tra1 Real Input tracors |
---|
80 | ! p1 Real Input full level pressure |
---|
81 | ! ph1 Real Input half level pressure |
---|
82 | ! ALE1 Real Input Available lifting Energy |
---|
83 | ! ALP1 Real Input Available lifting Power |
---|
84 | ! sig1feed1 Real Input sigma coord at lower bound of feeding layer |
---|
85 | ! sig2feed1 Real Input sigma coord at upper bound of feeding layer |
---|
86 | ! wght1 Real Input weight density determining the feeding mixture |
---|
87 | ! iflag1 Integer Output flag for Emanuel conditions |
---|
88 | ! ft1 Real Output temp tend |
---|
89 | ! fq1 Real Output spec hum tend |
---|
90 | ! fu1 Real Output u-wind tend |
---|
91 | ! fv1 Real Output v-wind tend |
---|
92 | ! ftra1 Real Output tracor tend |
---|
93 | ! precip1 Real Output precipitation |
---|
94 | ! kbas1 Integer Output cloud base level |
---|
95 | ! ktop1 Integer Output cloud top level |
---|
96 | ! cbmf1 Real Output cloud base mass flux |
---|
97 | ! sig1 Real In/Out section adiabatic updraft |
---|
98 | ! w01 Real In/Out vertical velocity within adiab updraft |
---|
99 | ! ptop21 Real In/Out top of entraining zone |
---|
100 | ! Ma1 Real Output mass flux adiabatic updraft |
---|
101 | ! mip1 Real Output mass flux shed by the adiabatic updraft |
---|
102 | ! Vprecip1 Real Output vertical profile of total precipitation |
---|
103 | ! Vprecipi1 Real Output vertical profile of ice precipitation |
---|
104 | ! upwd1 Real Output total upward mass flux (adiab+mixed) |
---|
105 | ! dnwd1 Real Output saturated downward mass flux (mixed) |
---|
106 | ! dnwd01 Real Output unsaturated downward mass flux |
---|
107 | ! qcondc1 Real Output in-cld mixing ratio of condensed water |
---|
108 | ! wd1 Real Output downdraft velocity scale for sfc fluxes |
---|
109 | ! cape1 Real Output CAPE |
---|
110 | ! cin1 Real Output CIN |
---|
111 | ! tvp1 Real Output adiab lifted parcell virt temp |
---|
112 | ! ftd1 Real Output precip temp tend |
---|
113 | ! fqt1 Real Output precip spec hum tend |
---|
114 | ! Plim11 Real Output |
---|
115 | ! Plim21 Real Output |
---|
116 | ! asupmax1 Real Output |
---|
117 | ! supmax01 Real Output |
---|
118 | ! asupmaxmin1 Real Output |
---|
119 | |
---|
120 | ! ftd1 Real Output Array of temperature tendency due to precipitations (K/s) of dimension ND, |
---|
121 | ! defined at same grid levels as T, Q, QS and P. |
---|
122 | |
---|
123 | ! fqd1 Real Output Array of specific humidity tendencies due to precipitations ((gm/gm)/s) |
---|
124 | ! of dimension ND, defined at same grid levels as T, Q, QS and P. |
---|
125 | |
---|
126 | ! wdtrainA1 Real Output precipitation detrained from adiabatic draught; |
---|
127 | ! used in tracer transport (cvltr) |
---|
128 | ! wdtrainM1 Real Output precipitation detrained from mixed draughts; |
---|
129 | ! used in tracer transport (cvltr) |
---|
130 | ! da1 Real Output used in tracer transport (cvltr) |
---|
131 | ! phi1 Real Output used in tracer transport (cvltr) |
---|
132 | ! mp1 Real Output used in tracer transport (cvltr) |
---|
133 | ! qtc1 Real Output specific humidity in convection |
---|
134 | ! sigt1 Real Output surface fraction in adiabatic updrafts |
---|
135 | ! phi21 Real Output used in tracer transport (cvltr) |
---|
136 | |
---|
137 | ! d1a1 Real Output used in tracer transport (cvltr) |
---|
138 | ! dam1 Real Output used in tracer transport (cvltr) |
---|
139 | |
---|
140 | ! epmlmMm1 Real Output used in tracer transport (cvltr) |
---|
141 | ! eplaMm1 Real Output used in tracer transport (cvltr) |
---|
142 | |
---|
143 | ! evap1 Real Output |
---|
144 | ! ep1 Real Output |
---|
145 | ! sigij1 Real Output used in tracer transport (cvltr) |
---|
146 | ! clw1 Real Output condensed water content of the adiabatic updraught |
---|
147 | ! elij1 Real Output |
---|
148 | ! wghti1 Real Output final weight of the feeding layers, |
---|
149 | ! used in tracer transport (cvltr) |
---|
150 | |
---|
151 | |
---|
152 | ! S. Bony, Mar 2002: |
---|
153 | ! * Several modules corresponding to different physical processes |
---|
154 | ! * Several versions of convect may be used: |
---|
155 | ! - iflag_con=3: version lmd (previously named convect3) |
---|
156 | ! - iflag_con=4: version 4.3b (vect. version, previously convect1/2) |
---|
157 | ! + tard: - iflag_con=5: version lmd with ice (previously named convectg) |
---|
158 | ! S. Bony, Oct 2002: |
---|
159 | ! * Vectorization of convect3 (ie version lmd) |
---|
160 | |
---|
161 | ! ..............................END PROLOGUE............................. |
---|
162 | |
---|
163 | |
---|
164 | |
---|
165 | ! Input |
---|
166 | INTEGER, INTENT (IN) :: len |
---|
167 | INTEGER, INTENT (IN) :: nd |
---|
168 | INTEGER, INTENT (IN) :: ndp1 |
---|
169 | INTEGER, INTENT (IN) :: ntra |
---|
170 | INTEGER, INTENT(IN) :: nloc ! (nloc=len) pour l'instant |
---|
171 | INTEGER, INTENT (IN) :: k_upper |
---|
172 | INTEGER, INTENT (IN) :: iflag_con |
---|
173 | INTEGER, INTENT (IN) :: iflag_mix |
---|
174 | INTEGER, INTENT (IN) :: iflag_ice_thermo |
---|
175 | INTEGER, INTENT (IN) :: iflag_clos |
---|
176 | LOGICAL, INTENT (IN) :: ok_conserv_q |
---|
177 | REAL, INTENT (IN) :: tau_cld_cv |
---|
178 | REAL, INTENT (IN) :: coefw_cld_cv |
---|
179 | REAL, INTENT (IN) :: delt |
---|
180 | REAL, INTENT (IN) :: comp_threshold |
---|
181 | REAL, DIMENSION (len, nd), INTENT (IN) :: t1 |
---|
182 | REAL, DIMENSION (len, nd), INTENT (IN) :: q1 |
---|
183 | REAL, DIMENSION (len, nd), INTENT (IN) :: qs1 |
---|
184 | REAL, DIMENSION (len, nd), INTENT (IN) :: t1_wake |
---|
185 | REAL, DIMENSION (len, nd), INTENT (IN) :: q1_wake |
---|
186 | REAL, DIMENSION (len, nd), INTENT (IN) :: qs1_wake |
---|
187 | REAL, DIMENSION (len), INTENT (IN) :: s1_wake |
---|
188 | REAL, DIMENSION (len, nd), INTENT (IN) :: u1 |
---|
189 | REAL, DIMENSION (len, nd), INTENT (IN) :: v1 |
---|
190 | REAL, DIMENSION (len, nd, ntra), INTENT (IN) :: tra1 |
---|
191 | REAL, DIMENSION (len, nd), INTENT (IN) :: p1 |
---|
192 | REAL, DIMENSION (len, ndp1), INTENT (IN) :: ph1 |
---|
193 | REAL, DIMENSION (len), INTENT (IN) :: Ale1 |
---|
194 | REAL, DIMENSION (len), INTENT (IN) :: Alp1 |
---|
195 | REAL, DIMENSION (len, nd), INTENT (IN) :: omega1 |
---|
196 | REAL, INTENT (IN) :: sig1feed1 ! pressure at lower bound of feeding layer |
---|
197 | REAL, INTENT (IN) :: sig2feed1 ! pressure at upper bound of feeding layer |
---|
198 | REAL, DIMENSION (nd), INTENT (IN) :: wght1 ! weight density determining the feeding mixture |
---|
199 | INTEGER, DIMENSION (len), INTENT (IN) :: lalim_conv1 |
---|
200 | |
---|
201 | ! Input/Output |
---|
202 | REAL, DIMENSION (len, nd), INTENT (INOUT) :: sig1 |
---|
203 | REAL, DIMENSION (len, nd), INTENT (INOUT) :: w01 |
---|
204 | |
---|
205 | ! Output |
---|
206 | INTEGER, DIMENSION (len), INTENT (OUT) :: iflag1 |
---|
207 | REAL, DIMENSION (len, nd), INTENT (OUT) :: ft1 |
---|
208 | REAL, DIMENSION (len, nd), INTENT (OUT) :: fq1 |
---|
209 | REAL, DIMENSION (len, nd), INTENT (OUT) :: fu1 |
---|
210 | REAL, DIMENSION (len, nd), INTENT (OUT) :: fv1 |
---|
211 | REAL, DIMENSION (len, nd, ntra), INTENT (OUT) :: ftra1 |
---|
212 | REAL, DIMENSION (len), INTENT (OUT) :: precip1 |
---|
213 | INTEGER, DIMENSION (len), INTENT (OUT) :: kbas1 |
---|
214 | INTEGER, DIMENSION (len), INTENT (OUT) :: ktop1 |
---|
215 | REAL, DIMENSION (len), INTENT (OUT) :: cbmf1 |
---|
216 | REAL, DIMENSION (len), INTENT (OUT) :: plcl1 |
---|
217 | REAL, DIMENSION (len), INTENT (OUT) :: plfc1 |
---|
218 | REAL, DIMENSION (len), INTENT (OUT) :: wbeff1 |
---|
219 | REAL, DIMENSION (len), INTENT (OUT) :: ptop21 |
---|
220 | REAL, DIMENSION (len), INTENT (OUT) :: sigd1 |
---|
221 | REAL, DIMENSION (len, nd), INTENT (OUT) :: ma1 ! adiab. asc. mass flux (staggered grid) |
---|
222 | REAL, DIMENSION (len, nd), INTENT (OUT) :: mip1 ! mass flux shed from adiab. ascent (extensive) |
---|
223 | ! real Vprecip1(len,nd) |
---|
224 | REAL, DIMENSION (len, ndp1), INTENT (OUT) :: vprecip1 ! tot precipitation flux (staggered grid) |
---|
225 | REAL, DIMENSION (len, ndp1), INTENT (OUT) :: vprecipi1 ! ice precipitation flux (staggered grid) |
---|
226 | REAL, DIMENSION (len, nd), INTENT (OUT) :: upwd1 ! upwd sat. mass flux (staggered grid) |
---|
227 | REAL, DIMENSION (len, nd), INTENT (OUT) :: dnwd1 ! dnwd sat. mass flux (staggered grid) |
---|
228 | REAL, DIMENSION (len, nd), INTENT (OUT) :: dnwd01 ! unsat. mass flux (staggered grid) |
---|
229 | REAL, DIMENSION (len, nd), INTENT (OUT) :: qcondc1 ! max cloud condensate (intensive) ! cld |
---|
230 | REAL, DIMENSION (len), INTENT (OUT) :: wd1 ! gust |
---|
231 | REAL, DIMENSION (len), INTENT (OUT) :: cape1 |
---|
232 | REAL, DIMENSION (len), INTENT (OUT) :: cin1 |
---|
233 | REAL, DIMENSION (len, nd), INTENT (OUT) :: tvp1 ! Virt. temp. in the adiab. ascent |
---|
234 | |
---|
235 | !AC! |
---|
236 | !! real da1(len,nd),phi1(len,nd,nd) |
---|
237 | !! real da(len,nd),phi(len,nd,nd) |
---|
238 | !AC! |
---|
239 | REAL, DIMENSION (len, nd), INTENT (OUT) :: ftd1 ! Temp. tendency due to the sole unsat. drafts |
---|
240 | REAL, DIMENSION (len, nd), INTENT (OUT) :: fqd1 ! Moist. tendency due to the sole unsat. drafts |
---|
241 | REAL, DIMENSION (len), INTENT (OUT) :: Plim11 |
---|
242 | REAL, DIMENSION (len), INTENT (OUT) :: Plim21 |
---|
243 | REAL, DIMENSION (len, nd), INTENT (OUT) :: asupmax1 ! Highest mixing fraction of mixed updraughts |
---|
244 | REAL, DIMENSION (len), INTENT (OUT) :: supmax01 |
---|
245 | REAL, DIMENSION (len), INTENT (OUT) :: asupmaxmin1 |
---|
246 | REAL, DIMENSION (len, nd), INTENT (OUT) :: qtc1 ! in cloud water content (intensive) ! cld |
---|
247 | REAL, DIMENSION (len, nd), INTENT (OUT) :: sigt1 ! fract. cloud area (intensive) ! cld |
---|
248 | |
---|
249 | ! RomP >>> |
---|
250 | REAL, DIMENSION (len, nd), INTENT (OUT) :: wdtrainA1, wdtrainM1 ! precipitation sources (extensive) |
---|
251 | REAL, DIMENSION (len, nd), INTENT (OUT) :: mp1 ! unsat. mass flux (staggered grid) |
---|
252 | REAL, DIMENSION (len, nd), INTENT (OUT) :: da1 ! detrained mass flux of adiab. asc. air (extensive) |
---|
253 | REAL, DIMENSION (len, nd, nd), INTENT (OUT) :: phi1 ! mass flux of envt. air in mixed draughts (extensive) |
---|
254 | REAL, DIMENSION (len, nd, nd), INTENT (OUT) :: epmlmMm1 ! (extensive) |
---|
255 | REAL, DIMENSION (len, nd), INTENT (OUT) :: eplaMm1 ! (extensive) |
---|
256 | REAL, DIMENSION (len, nd), INTENT (OUT) :: evap1 ! evaporation rate in precip. downdraft. (intensive) |
---|
257 | REAL, DIMENSION (len, nd), INTENT (OUT) :: ep1 |
---|
258 | REAL, DIMENSION (len, nd, nd), INTENT (OUT) :: sigij1 ! mass fraction of env. air in mixed draughts (intensive) |
---|
259 | REAL, DIMENSION (len, nd, nd), INTENT (OUT) :: elij1! cond. water per unit mass of mixed draughts (intensive) |
---|
260 | REAL, DIMENSION (len, nd), INTENT (OUT) :: clw1 ! cond. water per unit mass of the adiab. asc. (intensive) |
---|
261 | !JYG,RL |
---|
262 | REAL, DIMENSION (len, nd), INTENT (OUT) :: wghti1 ! final weight of the feeding layers (extensive) |
---|
263 | !JYG,RL |
---|
264 | REAL, DIMENSION (len, nd, nd), INTENT (OUT) :: phi21 ! (extensive) |
---|
265 | REAL, DIMENSION (len, nd), INTENT (OUT) :: d1a1 ! (extensive) |
---|
266 | REAL, DIMENSION (len, nd), INTENT (OUT) :: dam1 ! (extensive) |
---|
267 | ! RomP <<< |
---|
268 | REAL, DIMENSION (len ), INTENT (OUT) :: epmax_diag1 |
---|
269 | |
---|
270 | ! ------------------------------------------------------------------- |
---|
271 | ! Prolog by Kerry Emanuel. |
---|
272 | ! ------------------------------------------------------------------- |
---|
273 | ! --- ARGUMENTS |
---|
274 | ! ------------------------------------------------------------------- |
---|
275 | ! --- On input: |
---|
276 | |
---|
277 | ! t: Array of absolute temperature (K) of dimension ND, with first |
---|
278 | ! index corresponding to lowest model level. Note that this array |
---|
279 | ! will be altered by the subroutine if dry convective adjustment |
---|
280 | ! occurs and if IPBL is not equal to 0. |
---|
281 | |
---|
282 | ! q: Array of specific humidity (gm/gm) of dimension ND, with first |
---|
283 | ! index corresponding to lowest model level. Must be defined |
---|
284 | ! at same grid levels as T. Note that this array will be altered |
---|
285 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
286 | |
---|
287 | ! qs: Array of saturation specific humidity of dimension ND, with first |
---|
288 | ! index corresponding to lowest model level. Must be defined |
---|
289 | ! at same grid levels as T. Note that this array will be altered |
---|
290 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
291 | |
---|
292 | ! t_wake: Array of absolute temperature (K), seen by unsaturated draughts, |
---|
293 | ! of dimension ND, with first index corresponding to lowest model level. |
---|
294 | |
---|
295 | ! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts, |
---|
296 | ! of dimension ND, with first index corresponding to lowest model level. |
---|
297 | ! Must be defined at same grid levels as T. |
---|
298 | |
---|
299 | ! qs_wake: Array of saturation specific humidity, seen by unsaturated draughts, |
---|
300 | ! of dimension ND, with first index corresponding to lowest model level. |
---|
301 | ! Must be defined at same grid levels as T. |
---|
302 | |
---|
303 | ! s_wake: Array of fractionnal area occupied by the wakes. |
---|
304 | |
---|
305 | ! u: Array of zonal wind velocity (m/s) of dimension ND, witth first |
---|
306 | ! index corresponding with the lowest model level. Defined at |
---|
307 | ! same levels as T. Note that this array will be altered if |
---|
308 | ! dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
309 | |
---|
310 | ! v: Same as u but for meridional velocity. |
---|
311 | |
---|
312 | ! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA), |
---|
313 | ! where NTRA is the number of different tracers. If no |
---|
314 | ! convective tracer transport is needed, define a dummy |
---|
315 | ! input array of dimension (ND,1). Tracers are defined at |
---|
316 | ! same vertical levels as T. Note that this array will be altered |
---|
317 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
318 | |
---|
319 | ! p: Array of pressure (mb) of dimension ND, with first |
---|
320 | ! index corresponding to lowest model level. Must be defined |
---|
321 | ! at same grid levels as T. |
---|
322 | |
---|
323 | ! ph: Array of pressure (mb) of dimension ND+1, with first index |
---|
324 | ! corresponding to lowest level. These pressures are defined at |
---|
325 | ! levels intermediate between those of P, T, Q and QS. The first |
---|
326 | ! value of PH should be greater than (i.e. at a lower level than) |
---|
327 | ! the first value of the array P. |
---|
328 | |
---|
329 | ! ALE: Available lifting Energy |
---|
330 | |
---|
331 | ! ALP: Available lifting Power |
---|
332 | |
---|
333 | ! nl: The maximum number of levels to which convection can penetrate, plus 1. |
---|
334 | ! NL MUST be less than or equal to ND-1. |
---|
335 | |
---|
336 | ! delt: The model time step (sec) between calls to CONVECT |
---|
337 | |
---|
338 | ! ---------------------------------------------------------------------------- |
---|
339 | ! --- On Output: |
---|
340 | |
---|
341 | ! iflag: An output integer whose value denotes the following: |
---|
342 | ! VALUE INTERPRETATION |
---|
343 | ! ----- -------------- |
---|
344 | ! 0 Moist convection occurs. |
---|
345 | ! 1 Moist convection occurs, but a CFL condition |
---|
346 | ! on the subsidence warming is violated. This |
---|
347 | ! does not cause the scheme to terminate. |
---|
348 | ! 2 Moist convection, but no precip because ep(inb) lt 0.0001 |
---|
349 | ! 3 No moist convection because new cbmf is 0 and old cbmf is 0. |
---|
350 | ! 4 No moist convection; atmosphere is not |
---|
351 | ! unstable |
---|
352 | ! 6 No moist convection because ihmin le minorig. |
---|
353 | ! 7 No moist convection because unreasonable |
---|
354 | ! parcel level temperature or specific humidity. |
---|
355 | ! 8 No moist convection: lifted condensation |
---|
356 | ! level is above the 200 mb level. |
---|
357 | ! 9 No moist convection: cloud base is higher |
---|
358 | ! then the level NL-1. |
---|
359 | ! 10 No moist convection: cloud top is too warm. |
---|
360 | ! |
---|
361 | |
---|
362 | ! ft: Array of temperature tendency (K/s) of dimension ND, defined at same |
---|
363 | ! grid levels as T, Q, QS and P. |
---|
364 | |
---|
365 | ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND, |
---|
366 | ! defined at same grid levels as T, Q, QS and P. |
---|
367 | |
---|
368 | ! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND, |
---|
369 | ! defined at same grid levels as T. |
---|
370 | |
---|
371 | ! fv: Same as FU, but for forcing of meridional velocity. |
---|
372 | |
---|
373 | ! ftra: Array of forcing of tracer content, in tracer mixing ratio per |
---|
374 | ! second, defined at same levels as T. Dimensioned (ND,NTRA). |
---|
375 | |
---|
376 | ! precip: Scalar convective precipitation rate (mm/day). |
---|
377 | |
---|
378 | ! wd: A convective downdraft velocity scale. For use in surface |
---|
379 | ! flux parameterizations. See convect.ps file for details. |
---|
380 | |
---|
381 | ! tprime: A convective downdraft temperature perturbation scale (K). |
---|
382 | ! For use in surface flux parameterizations. See convect.ps |
---|
383 | ! file for details. |
---|
384 | |
---|
385 | ! qprime: A convective downdraft specific humidity |
---|
386 | ! perturbation scale (gm/gm). |
---|
387 | ! For use in surface flux parameterizations. See convect.ps |
---|
388 | ! file for details. |
---|
389 | |
---|
390 | ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST |
---|
391 | ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT |
---|
392 | ! ITS NEXT CALL. That is, the value of CBMF must be "remembered" |
---|
393 | ! by the calling program between calls to CONVECT. |
---|
394 | |
---|
395 | ! det: Array of detrainment mass flux of dimension ND. |
---|
396 | ! ------------------------------------------------------------------- |
---|
397 | |
---|
398 | ! Local (non compressed) arrays |
---|
399 | |
---|
400 | |
---|
401 | INTEGER i, k, il |
---|
402 | INTEGER nword1, nword2, nword3, nword4 |
---|
403 | INTEGER icbmax |
---|
404 | INTEGER nk1(len) |
---|
405 | INTEGER icb1(len) |
---|
406 | INTEGER icbs1(len) |
---|
407 | |
---|
408 | LOGICAL ok_inhib ! True => possible inhibition of convection by dryness |
---|
409 | LOGICAL, SAVE :: debut = .TRUE. |
---|
410 | !$OMP THREADPRIVATE(debut) |
---|
411 | |
---|
412 | REAL coef_convective(len) ! = 1 for convective points, = 0 otherwise |
---|
413 | REAL tnk1(len) |
---|
414 | REAL thnk1(len) |
---|
415 | REAL qnk1(len) |
---|
416 | REAL gznk1(len) |
---|
417 | REAL qsnk1(len) |
---|
418 | REAL unk1(len) |
---|
419 | REAL vnk1(len) |
---|
420 | REAL cpnk1(len) |
---|
421 | REAL hnk1(len) |
---|
422 | REAL pbase1(len) |
---|
423 | REAL buoybase1(len) |
---|
424 | |
---|
425 | REAL lf1(len, nd), lf1_wake(len, nd) |
---|
426 | REAL lv1(len, nd), lv1_wake(len, nd) |
---|
427 | REAL cpn1(len, nd), cpn1_wake(len, nd) |
---|
428 | REAL tv1(len, nd), tv1_wake(len, nd) |
---|
429 | REAL gz1(len, nd), gz1_wake(len, nd) |
---|
430 | REAL hm1(len, nd) |
---|
431 | REAL h1(len, nd), h1_wake(len, nd) |
---|
432 | REAL tp1(len, nd) |
---|
433 | REAL th1(len, nd), th1_wake(len, nd) |
---|
434 | |
---|
435 | REAL bid(len, nd) ! dummy array |
---|
436 | |
---|
437 | INTEGER ncum |
---|
438 | |
---|
439 | REAL p1feed1(len) ! pressure at lower bound of feeding layer |
---|
440 | REAL p2feed1(len) ! pressure at upper bound of feeding layer |
---|
441 | !JYG,RL |
---|
442 | !! real wghti1(len,nd) ! weights of the feeding layers |
---|
443 | !JYG,RL |
---|
444 | |
---|
445 | ! (local) compressed fields: |
---|
446 | |
---|
447 | |
---|
448 | INTEGER idcum(nloc) |
---|
449 | !jyg< |
---|
450 | LOGICAL compress ! True if compression occurs |
---|
451 | !>jyg |
---|
452 | INTEGER iflag(nloc), nk(nloc), icb(nloc) |
---|
453 | INTEGER nent(nloc, nd) |
---|
454 | INTEGER icbs(nloc) |
---|
455 | INTEGER inb(nloc), inbis(nloc) |
---|
456 | |
---|
457 | REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc) |
---|
458 | REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd) |
---|
459 | REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd) |
---|
460 | REAL s_wake(nloc) |
---|
461 | REAL u(nloc, nd), v(nloc, nd) |
---|
462 | REAL gz(nloc, nd), h(nloc, nd) |
---|
463 | REAL h_wake(nloc, nd) |
---|
464 | REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd) |
---|
465 | REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd) |
---|
466 | REAL p(nloc, nd), ph(nloc, nd+1), tv(nloc, nd), tp(nloc, nd) |
---|
467 | REAL tv_wake(nloc, nd) |
---|
468 | REAL clw(nloc, nd) |
---|
469 | REAL dph(nloc, nd) |
---|
470 | REAL pbase(nloc), buoybase(nloc), th(nloc, nd) |
---|
471 | REAL th_wake(nloc, nd) |
---|
472 | REAL tvp(nloc, nd) |
---|
473 | REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc) |
---|
474 | REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd) |
---|
475 | REAL buoy(nloc, nd) |
---|
476 | REAL cape(nloc) |
---|
477 | REAL cin(nloc) |
---|
478 | REAL m(nloc, nd) |
---|
479 | REAL ment(nloc, nd, nd), sigij(nloc, nd, nd) |
---|
480 | REAL qent(nloc, nd, nd) |
---|
481 | REAL hent(nloc, nd, nd) |
---|
482 | REAL uent(nloc, nd, nd), vent(nloc, nd, nd) |
---|
483 | REAL ments(nloc, nd, nd), qents(nloc, nd, nd) |
---|
484 | REAL elij(nloc, nd, nd) |
---|
485 | REAL supmax(nloc, nd) |
---|
486 | REAL Ale(nloc), Alp(nloc), coef_clos(nloc) |
---|
487 | REAL omega(nloc,nd) |
---|
488 | REAL sigd(nloc) |
---|
489 | ! real mp(nloc,nd), qp(nloc,nd), up(nloc,nd), vp(nloc,nd) |
---|
490 | ! real wt(nloc,nd), water(nloc,nd), evap(nloc,nd), ice(nloc,nd) |
---|
491 | ! real b(nloc,nd), sigd(nloc) |
---|
492 | ! save mp,qp,up,vp,wt,water,evap,b |
---|
493 | REAL, DIMENSION(len,nd) :: mp, qp, up, vp |
---|
494 | REAL, DIMENSION(len,nd) :: wt, water, evap |
---|
495 | REAL, DIMENSION(len,nd) :: ice, fondue, b |
---|
496 | REAL, DIMENSION(len,nd) :: frac, faci |
---|
497 | REAL ft(nloc, nd), fq(nloc, nd) |
---|
498 | REAL ftd(nloc, nd), fqd(nloc, nd) |
---|
499 | REAL fu(nloc, nd), fv(nloc, nd) |
---|
500 | REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd) |
---|
501 | REAL ma(nloc, nd), mip(nloc, nd) |
---|
502 | !! REAL tls(nloc, nd), tps(nloc, nd) ! unused . jyg |
---|
503 | REAL qprime(nloc), tprime(nloc) |
---|
504 | REAL precip(nloc) |
---|
505 | ! real Vprecip(nloc,nd) |
---|
506 | REAL vprecip(nloc, nd+1) |
---|
507 | REAL vprecipi(nloc, nd+1) |
---|
508 | REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra) |
---|
509 | REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra) |
---|
510 | REAL qcondc(nloc, nd) ! cld |
---|
511 | REAL wd(nloc) ! gust |
---|
512 | REAL Plim1(nloc), plim2(nloc) |
---|
513 | REAL asupmax(nloc, nd) |
---|
514 | REAL supmax0(nloc) |
---|
515 | REAL asupmaxmin(nloc) |
---|
516 | |
---|
517 | REAL tnk(nloc), qnk(nloc), gznk(nloc) |
---|
518 | REAL wghti(nloc, nd) |
---|
519 | REAL hnk(nloc), unk(nloc), vnk(nloc) |
---|
520 | |
---|
521 | REAL qtc(nloc, nd) ! cld |
---|
522 | REAL sigt(nloc, nd) ! cld |
---|
523 | |
---|
524 | ! RomP >>> |
---|
525 | REAL wdtrainA(nloc, nd), wdtrainM(nloc, nd) |
---|
526 | REAL da(len, nd), phi(len, nd, nd) |
---|
527 | REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd) |
---|
528 | REAL phi2(len, nd, nd) |
---|
529 | REAL d1a(len, nd), dam(len, nd) |
---|
530 | ! RomP <<< |
---|
531 | REAL epmax_diag(nloc) ! epmax_cape |
---|
532 | |
---|
533 | CHARACTER (LEN=20) :: modname = 'cva_driver' |
---|
534 | CHARACTER (LEN=80) :: abort_message |
---|
535 | |
---|
536 | INTEGER,SAVE :: igout=1 |
---|
537 | !$OMP THREADPRIVATE(igout) |
---|
538 | |
---|
539 | |
---|
540 | ! print *, 't1, t1_wake ',(k,t1(1,k),t1_wake(1,k),k=1,nd) |
---|
541 | ! print *, 'q1, q1_wake ',(k,q1(1,k),q1_wake(1,k),k=1,nd) |
---|
542 | |
---|
543 | ! ------------------------------------------------------------------- |
---|
544 | ! --- SET CONSTANTS AND PARAMETERS |
---|
545 | ! ------------------------------------------------------------------- |
---|
546 | |
---|
547 | ! -- set simulation flags: |
---|
548 | ! (common cvflag) |
---|
549 | |
---|
550 | CALL cv_flag(iflag_ice_thermo) |
---|
551 | |
---|
552 | ! -- set thermodynamical constants: |
---|
553 | ! (common cvthermo) |
---|
554 | |
---|
555 | CALL cv_thermo(iflag_con) |
---|
556 | |
---|
557 | ! -- set convect parameters |
---|
558 | |
---|
559 | ! includes microphysical parameters and parameters that |
---|
560 | ! control the rate of approach to quasi-equilibrium) |
---|
561 | ! (common cvparam) |
---|
562 | |
---|
563 | IF (iflag_con==3) THEN |
---|
564 | CALL cv3_param(nd, k_upper, delt) |
---|
565 | |
---|
566 | END IF |
---|
567 | |
---|
568 | IF (iflag_con==4) THEN |
---|
569 | CALL cv_param(nd) |
---|
570 | END IF |
---|
571 | |
---|
572 | ! --------------------------------------------------------------------- |
---|
573 | ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS |
---|
574 | ! --------------------------------------------------------------------- |
---|
575 | nword1 = len |
---|
576 | nword2 = len*nd |
---|
577 | nword3 = len*nd*ntra |
---|
578 | nword4 = len*nd*nd |
---|
579 | |
---|
580 | iflag1(:) = 0 |
---|
581 | ktop1(:) = 0 |
---|
582 | kbas1(:) = 0 |
---|
583 | ft1(:, :) = 0.0 |
---|
584 | fq1(:, :) = 0.0 |
---|
585 | fu1(:, :) = 0.0 |
---|
586 | fv1(:, :) = 0.0 |
---|
587 | ftra1(:, :, :) = 0. |
---|
588 | precip1(:) = 0. |
---|
589 | cbmf1(:) = 0. |
---|
590 | plcl1(:) = 0. |
---|
591 | plfc1(:) = 0. |
---|
592 | wbeff1(:) = 0. |
---|
593 | ptop21(:) = 0. |
---|
594 | sigd1(:) = 0. |
---|
595 | ma1(:, :) = 0. |
---|
596 | mip1(:, :) = 0. |
---|
597 | vprecip1(:, :) = 0. |
---|
598 | vprecipi1(:, :) = 0. |
---|
599 | upwd1(:, :) = 0. |
---|
600 | dnwd1(:, :) = 0. |
---|
601 | dnwd01(:, :) = 0. |
---|
602 | qcondc1(:, :) = 0. |
---|
603 | wd1(:) = 0. |
---|
604 | cape1(:) = 0. |
---|
605 | cin1(:) = 0. |
---|
606 | tvp1(:, :) = 0. |
---|
607 | ftd1(:, :) = 0. |
---|
608 | fqd1(:, :) = 0. |
---|
609 | Plim11(:) = 0. |
---|
610 | Plim21(:) = 0. |
---|
611 | asupmax1(:, :) = 0. |
---|
612 | supmax01(:) = 0. |
---|
613 | asupmaxmin1(:) = 0. |
---|
614 | |
---|
615 | DO il = 1, len |
---|
616 | cin1(il) = -100000. |
---|
617 | cape1(il) = -1. |
---|
618 | END DO |
---|
619 | |
---|
620 | !! IF (iflag_con==3) THEN |
---|
621 | !! DO il = 1, len |
---|
622 | !! sig1(il, nd) = sig1(il, nd) + 1. |
---|
623 | !! sig1(il, nd) = amin1(sig1(il,nd), 12.1) |
---|
624 | !! END DO |
---|
625 | !! END IF |
---|
626 | |
---|
627 | IF (iflag_con==3) THEN |
---|
628 | CALL cv3_incrcount(len,nd,delt,sig1) |
---|
629 | END IF ! (iflag_con==3) |
---|
630 | |
---|
631 | ! RomP >>> |
---|
632 | sigt1(:, :) = 0. |
---|
633 | qtc1(:, :) = 0. |
---|
634 | wdtrainA1(:, :) = 0. |
---|
635 | wdtrainM1(:, :) = 0. |
---|
636 | da1(:, :) = 0. |
---|
637 | phi1(:, :, :) = 0. |
---|
638 | epmlmMm1(:, :, :) = 0. |
---|
639 | eplaMm1(:, :) = 0. |
---|
640 | mp1(:, :) = 0. |
---|
641 | evap1(:, :) = 0. |
---|
642 | ep1(:, :) = 0. |
---|
643 | sigij1(:, :, :) = 0. |
---|
644 | elij1(:, :, :) = 0. |
---|
645 | clw1(:,:) = 0. |
---|
646 | wghti1(:,:) = 0. |
---|
647 | phi21(:, :, :) = 0. |
---|
648 | d1a1(:, :) = 0. |
---|
649 | dam1(:, :) = 0. |
---|
650 | ! RomP <<< |
---|
651 | ! --------------------------------------------------------------------- |
---|
652 | ! --- INITIALIZE LOCAL ARRAYS AND PARAMETERS |
---|
653 | ! --------------------------------------------------------------------- |
---|
654 | |
---|
655 | DO il = 1, nloc |
---|
656 | coef_clos(il) = 1. |
---|
657 | END DO |
---|
658 | |
---|
659 | ! -------------------------------------------------------------------- |
---|
660 | ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY |
---|
661 | ! -------------------------------------------------------------------- |
---|
662 | |
---|
663 | IF (iflag_con==3) THEN |
---|
664 | |
---|
665 | IF (debut) THEN |
---|
666 | PRINT *, 'Emanuel version 3 nouvelle' |
---|
667 | END IF |
---|
668 | ! print*,'t1, q1 ',t1,q1 |
---|
669 | if (prt_level >= 9) & |
---|
670 | PRINT *, 'cva_driver -> cv3_prelim' |
---|
671 | CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, & ! nd->na |
---|
672 | lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1) |
---|
673 | |
---|
674 | |
---|
675 | if (prt_level >= 9) & |
---|
676 | PRINT *, 'cva_driver -> cv3_prelim' |
---|
677 | CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, & ! nd->na |
---|
678 | lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, & |
---|
679 | h1_wake, bid, th1_wake) |
---|
680 | |
---|
681 | END IF |
---|
682 | |
---|
683 | IF (iflag_con==4) THEN |
---|
684 | PRINT *, 'Emanuel version 4 ' |
---|
685 | if (prt_level >= 9) & |
---|
686 | PRINT *, 'cva_driver -> cv_prelim' |
---|
687 | CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1, & |
---|
688 | lv1, cpn1, tv1, gz1, h1, hm1) |
---|
689 | END IF |
---|
690 | |
---|
691 | ! -------------------------------------------------------------------- |
---|
692 | ! --- CONVECTIVE FEED |
---|
693 | ! -------------------------------------------------------------------- |
---|
694 | |
---|
695 | ! compute feeding layer potential temperature and mixing ratio : |
---|
696 | |
---|
697 | ! get bounds of feeding layer |
---|
698 | |
---|
699 | ! test niveaux couche alimentation KE |
---|
700 | IF (sig1feed1==sig2feed1) THEN |
---|
701 | WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed' |
---|
702 | WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def' |
---|
703 | abort_message = '' |
---|
704 | CALL abort_physic(modname, abort_message, 1) |
---|
705 | END IF |
---|
706 | |
---|
707 | DO i = 1, len |
---|
708 | p1feed1(i) = sig1feed1*ph1(i, 1) |
---|
709 | p2feed1(i) = sig2feed1*ph1(i, 1) |
---|
710 | !test maf |
---|
711 | ! p1feed1(i)=ph1(i,1) |
---|
712 | ! p2feed1(i)=ph1(i,2) |
---|
713 | ! p2feed1(i)=ph1(i,3) |
---|
714 | !testCR: on prend la couche alim des thermiques |
---|
715 | ! p2feed1(i)=ph1(i,lalim_conv1(i)+1) |
---|
716 | ! print*,'lentr=',lentr(i),ph1(i,lentr(i)+1),ph1(i,2) |
---|
717 | END DO |
---|
718 | |
---|
719 | IF (iflag_con==3) THEN |
---|
720 | END IF |
---|
721 | DO i = 1, len |
---|
722 | ! print*,'avant cv3_feed Plim',p1feed1(i),p2feed1(i) |
---|
723 | END DO |
---|
724 | IF (iflag_con==3) THEN |
---|
725 | |
---|
726 | ! print*, 'IFLAG1 avant cv3_feed' |
---|
727 | ! print*,'len,nd',len,nd |
---|
728 | ! write(*,'(64i1)') iflag1(2:len-1) |
---|
729 | |
---|
730 | if (prt_level >= 9) & |
---|
731 | PRINT *, 'cva_driver -> cv3_feed' |
---|
732 | CALL cv3_feed(len, nd, ok_conserv_q, & ! nd->na |
---|
733 | t1, q1, u1, v1, p1, ph1, h1, gz1, & |
---|
734 | p1feed1, p2feed1, wght1, & |
---|
735 | wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, & |
---|
736 | cpnk1, hnk1, nk1, icb1, icbmax, iflag1, gznk1, plcl1) |
---|
737 | END IF |
---|
738 | |
---|
739 | ! print*, 'IFLAG1 apres cv3_feed' |
---|
740 | ! print*,'len,nd',len,nd |
---|
741 | ! write(*,'(64i1)') iflag1(2:len-1) |
---|
742 | |
---|
743 | IF (iflag_con==4) THEN |
---|
744 | if (prt_level >= 9) & |
---|
745 | PRINT *, 'cva_driver -> cv_feed' |
---|
746 | CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1, & |
---|
747 | nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) |
---|
748 | END IF |
---|
749 | |
---|
750 | ! print *, 'cv3_feed-> iflag1, plcl1 ',iflag1(1),plcl1(1) |
---|
751 | |
---|
752 | ! -------------------------------------------------------------------- |
---|
753 | ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part |
---|
754 | ! (up through ICB for convect4, up through ICB+1 for convect3) |
---|
755 | ! Calculates the lifted parcel virtual temperature at nk, the |
---|
756 | ! actual temperature, and the adiabatic liquid water content. |
---|
757 | ! -------------------------------------------------------------------- |
---|
758 | |
---|
759 | IF (iflag_con==3) THEN |
---|
760 | |
---|
761 | if (prt_level >= 9) & |
---|
762 | PRINT *, 'cva_driver -> cv3_undilute1' |
---|
763 | CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, icb1, tnk1, qnk1, & ! nd->na |
---|
764 | gznk1, tp1, tvp1, clw1, icbs1) |
---|
765 | END IF |
---|
766 | |
---|
767 | |
---|
768 | IF (iflag_con==4) THEN |
---|
769 | if (prt_level >= 9) & |
---|
770 | PRINT *, 'cva_driver -> cv_undilute1' |
---|
771 | CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, & |
---|
772 | tp1, tvp1, clw1) |
---|
773 | END IF |
---|
774 | |
---|
775 | ! ------------------------------------------------------------------- |
---|
776 | ! --- TRIGGERING |
---|
777 | ! ------------------------------------------------------------------- |
---|
778 | |
---|
779 | ! print *,' avant triggering, iflag_con ',iflag_con |
---|
780 | |
---|
781 | IF (iflag_con==3) THEN |
---|
782 | |
---|
783 | if (prt_level >= 9) & |
---|
784 | PRINT *, 'cva_driver -> cv3_trigger' |
---|
785 | CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na |
---|
786 | pbase1, buoybase1, iflag1, sig1, w01) |
---|
787 | |
---|
788 | |
---|
789 | ! print*, 'IFLAG1 apres cv3_triger' |
---|
790 | ! print*,'len,nd',len,nd |
---|
791 | ! write(*,'(64i1)') iflag1(2:len-1) |
---|
792 | |
---|
793 | ! call dump2d(iim,jjm-1,sig1(2) |
---|
794 | END IF |
---|
795 | |
---|
796 | IF (iflag_con==4) THEN |
---|
797 | if (prt_level >= 9) & |
---|
798 | PRINT *, 'cva_driver -> cv_trigger' |
---|
799 | CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1) |
---|
800 | END IF |
---|
801 | |
---|
802 | |
---|
803 | ! ===================================================================== |
---|
804 | ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY |
---|
805 | ! ===================================================================== |
---|
806 | |
---|
807 | ! Determine the number "ncum" of convective gridpoints, the list "idcum" of convective |
---|
808 | ! gridpoints and the weights "coef_convective" (= 1. for convective gridpoints and = 0. |
---|
809 | ! elsewhere). |
---|
810 | ncum = 0 |
---|
811 | coef_convective(:) = 0. |
---|
812 | DO i = 1, len |
---|
813 | IF (iflag1(i)==0) THEN |
---|
814 | coef_convective(i) = 1. |
---|
815 | ncum = ncum + 1 |
---|
816 | idcum(ncum) = i |
---|
817 | END IF |
---|
818 | END DO |
---|
819 | |
---|
820 | ! print*,'len, ncum = ',len,ncum |
---|
821 | |
---|
822 | IF (ncum>0) THEN |
---|
823 | |
---|
824 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
825 | ! --- COMPRESS THE FIELDS |
---|
826 | ! (-> vectorization over convective gridpoints) |
---|
827 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
828 | |
---|
829 | IF (iflag_con==3) THEN |
---|
830 | ! print*,'ncum tv1 ',ncum,tv1 |
---|
831 | ! print*,'tvp1 ',tvp1 |
---|
832 | !jyg< |
---|
833 | ! If the fraction of convective points is larger than comp_threshold, then compression |
---|
834 | ! is assumed useless. |
---|
835 | ! |
---|
836 | compress = ncum .lt. len*comp_threshold |
---|
837 | ! |
---|
838 | IF (.not. compress) THEN |
---|
839 | DO i = 1,len |
---|
840 | idcum(i) = i |
---|
841 | ENDDO |
---|
842 | ENDIF |
---|
843 | ! |
---|
844 | !>jyg |
---|
845 | if (prt_level >= 9) & |
---|
846 | PRINT *, 'cva_driver -> cv3a_compress' |
---|
847 | CALL cv3a_compress(len, nloc, ncum, nd, ntra, compress, & |
---|
848 | iflag1, nk1, icb1, icbs1, & |
---|
849 | plcl1, tnk1, qnk1, gznk1, hnk1, unk1, vnk1, & |
---|
850 | wghti1, pbase1, buoybase1, & |
---|
851 | t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, & |
---|
852 | u1, v1, gz1, th1, th1_wake, & |
---|
853 | tra1, & |
---|
854 | h1, lv1, lf1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, & |
---|
855 | h1_wake, lv1_wake, lf1_wake, cpn1_wake, tv1_wake, & |
---|
856 | sig1, w01, ptop21, & |
---|
857 | Ale1, Alp1, omega1, & |
---|
858 | iflag, nk, icb, icbs, & |
---|
859 | plcl, tnk, qnk, gznk, hnk, unk, vnk, & |
---|
860 | wghti, pbase, buoybase, & |
---|
861 | t, q, qs, t_wake, q_wake, qs_wake, s_wake, & |
---|
862 | u, v, gz, th, th_wake, & |
---|
863 | tra, & |
---|
864 | h, lv, lf, cpn, p, ph, tv, tp, tvp, clw, & |
---|
865 | h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, & |
---|
866 | sig, w0, ptop2, & |
---|
867 | Ale, Alp, omega) |
---|
868 | |
---|
869 | ! print*,'tv ',tv |
---|
870 | ! print*,'tvp ',tvp |
---|
871 | |
---|
872 | END IF |
---|
873 | |
---|
874 | IF (iflag_con==4) THEN |
---|
875 | if (prt_level >= 9) & |
---|
876 | PRINT *, 'cva_driver -> cv_compress' |
---|
877 | CALL cv_compress(len, nloc, ncum, nd, & |
---|
878 | iflag1, nk1, icb1, & |
---|
879 | cbmf1, plcl1, tnk1, qnk1, gznk1, & |
---|
880 | t1, q1, qs1, u1, v1, gz1, & |
---|
881 | h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, & |
---|
882 | iflag, nk, icb, & |
---|
883 | cbmf, plcl, tnk, qnk, gznk, & |
---|
884 | t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw, & |
---|
885 | dph) |
---|
886 | END IF |
---|
887 | |
---|
888 | ! ------------------------------------------------------------------- |
---|
889 | ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part : |
---|
890 | ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
---|
891 | ! --- & |
---|
892 | ! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE |
---|
893 | ! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD |
---|
894 | ! --- & |
---|
895 | ! --- FIND THE LEVEL OF NEUTRAL BUOYANCY |
---|
896 | ! ------------------------------------------------------------------- |
---|
897 | |
---|
898 | IF (iflag_con==3) THEN |
---|
899 | if (prt_level >= 9) & |
---|
900 | PRINT *, 'cva_driver -> cv3_undilute2' |
---|
901 | CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, & !na->nd |
---|
902 | tnk, qnk, gznk, hnk, t, q, qs, gz, & |
---|
903 | p, ph, h, tv, lv, lf, pbase, buoybase, plcl, & |
---|
904 | inb, tp, tvp, clw, hp, ep, sigp, buoy, & |
---|
905 | frac) |
---|
906 | END IF |
---|
907 | |
---|
908 | IF (iflag_con==4) THEN |
---|
909 | if (prt_level >= 9) & |
---|
910 | PRINT *, 'cva_driver -> cv_undilute2' |
---|
911 | CALL cv_undilute2(nloc, ncum, nd, icb, nk, & |
---|
912 | tnk, qnk, gznk, t, q, qs, gz, & |
---|
913 | p, dph, h, tv, lv, & |
---|
914 | inb, inbis, tp, tvp, clw, hp, ep, sigp, frac) |
---|
915 | END IF |
---|
916 | |
---|
917 | ! epmax_cape |
---|
918 | ! on recalcule ep et hp |
---|
919 | if (prt_level >= 9) & |
---|
920 | PRINT *, 'cva_driver -> cv3_epmax_cape' |
---|
921 | call cv3_epmax_fn_cape(nloc,ncum,nd & |
---|
922 | , ep,hp,icb,inb,clw,nk,t,h,hnk,lv,lf,frac & |
---|
923 | , pbase, p, ph, tv, buoy, sig, w0,iflag & |
---|
924 | , epmax_diag) |
---|
925 | |
---|
926 | ! ------------------------------------------------------------------- |
---|
927 | ! --- MIXING(1) (if iflag_mix .ge. 1) |
---|
928 | ! ------------------------------------------------------------------- |
---|
929 | IF (iflag_con==3) THEN |
---|
930 | ! IF ((iflag_ice_thermo==1) .AND. (iflag_mix/=0)) THEN |
---|
931 | ! WRITE (*, *) ' iflag_ice_thermo==1 requires iflag_mix==0', ' but iflag_mix=', iflag_mix, & |
---|
932 | ! '. Might as well stop here.' |
---|
933 | ! STOP |
---|
934 | ! END IF |
---|
935 | IF (iflag_mix>=1) THEN |
---|
936 | CALL zilch(supmax, nloc*nd) |
---|
937 | if (prt_level >= 9) & |
---|
938 | PRINT *, 'cva_driver -> cv3p_mixing' |
---|
939 | CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, & ! na->nd |
---|
940 | ph, t, q, qs, u, v, tra, h, lv, lf, frac, qnk, & |
---|
941 | unk, vnk, hp, tv, tvp, ep, clw, sig, & |
---|
942 | ment, qent, hent, uent, vent, nent, & |
---|
943 | sigij, elij, supmax, ments, qents, traent) |
---|
944 | ! print*, 'cv3p_mixing-> supmax ', (supmax(1,k), k=1,nd) |
---|
945 | |
---|
946 | ELSE |
---|
947 | CALL zilch(supmax, nloc*nd) |
---|
948 | END IF |
---|
949 | END IF |
---|
950 | ! ------------------------------------------------------------------- |
---|
951 | ! --- CLOSURE |
---|
952 | ! ------------------------------------------------------------------- |
---|
953 | |
---|
954 | |
---|
955 | IF (iflag_con==3) THEN |
---|
956 | IF (iflag_clos==0) THEN |
---|
957 | if (prt_level >= 9) & |
---|
958 | PRINT *, 'cva_driver -> cv3_closure' |
---|
959 | CALL cv3_closure(nloc, ncum, nd, icb, inb, & ! na->nd |
---|
960 | pbase, p, ph, tv, buoy, & |
---|
961 | sig, w0, cape, m, iflag) |
---|
962 | END IF ! iflag_clos==0 |
---|
963 | |
---|
964 | ok_inhib = iflag_mix == 2 |
---|
965 | |
---|
966 | IF (iflag_clos==1) THEN |
---|
967 | PRINT *, ' pas d appel cv3p_closure' |
---|
968 | ! c CALL cv3p_closure(nloc,ncum,nd,icb,inb ! na->nd |
---|
969 | ! c : ,pbase,plcl,p,ph,tv,tvp,buoy |
---|
970 | ! c : ,supmax |
---|
971 | ! c o ,sig,w0,ptop2,cape,cin,m) |
---|
972 | END IF ! iflag_clos==1 |
---|
973 | |
---|
974 | IF (iflag_clos==2) THEN |
---|
975 | if (prt_level >= 9) & |
---|
976 | PRINT *, 'cva_driver -> cv3p1_closure' |
---|
977 | CALL cv3p1_closure(nloc, ncum, nd, icb, inb, & ! na->nd |
---|
978 | pbase, plcl, p, ph, tv, tvp, buoy, & |
---|
979 | supmax, ok_inhib, Ale, Alp, omega, & |
---|
980 | sig, w0, ptop2, cape, cin, m, iflag, coef_clos, & |
---|
981 | Plim1, plim2, asupmax, supmax0, & |
---|
982 | asupmaxmin, cbmf, plfc, wbeff) |
---|
983 | if (prt_level >= 10) & |
---|
984 | PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1) |
---|
985 | END IF ! iflag_clos==2 |
---|
986 | |
---|
987 | IF (iflag_clos==3) THEN |
---|
988 | if (prt_level >= 9) & |
---|
989 | PRINT *, 'cva_driver -> cv3p2_closure' |
---|
990 | CALL cv3p2_closure(nloc, ncum, nd, icb, inb, & ! na->nd |
---|
991 | pbase, plcl, p, ph, tv, tvp, buoy, & |
---|
992 | supmax, ok_inhib, Ale, Alp, omega, & |
---|
993 | sig, w0, ptop2, cape, cin, m, iflag, coef_clos, & |
---|
994 | Plim1, plim2, asupmax, supmax0, & |
---|
995 | asupmaxmin, cbmf, plfc, wbeff) |
---|
996 | if (prt_level >= 10) & |
---|
997 | PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1) |
---|
998 | END IF ! iflag_clos==3 |
---|
999 | END IF ! iflag_con==3 |
---|
1000 | |
---|
1001 | IF (iflag_con==4) THEN |
---|
1002 | if (prt_level >= 9) & |
---|
1003 | PRINT *, 'cva_driver -> cv_closure' |
---|
1004 | CALL cv_closure(nloc, ncum, nd, nk, icb, & |
---|
1005 | tv, tvp, p, ph, dph, plcl, cpn, & |
---|
1006 | iflag, cbmf) |
---|
1007 | END IF |
---|
1008 | |
---|
1009 | ! print *,'cv_closure-> cape ',cape(1) |
---|
1010 | |
---|
1011 | ! ------------------------------------------------------------------- |
---|
1012 | ! --- MIXING(2) |
---|
1013 | ! ------------------------------------------------------------------- |
---|
1014 | |
---|
1015 | IF (iflag_con==3) THEN |
---|
1016 | IF (iflag_mix==0) THEN |
---|
1017 | if (prt_level >= 9) & |
---|
1018 | PRINT *, 'cva_driver -> cv3_mixing' |
---|
1019 | CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, & ! na->nd |
---|
1020 | ph, t, q, qs, u, v, tra, h, lv, lf, frac, qnk, & |
---|
1021 | unk, vnk, hp, tv, tvp, ep, clw, m, sig, & |
---|
1022 | ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent) |
---|
1023 | CALL zilch(hent, nloc*nd*nd) |
---|
1024 | ELSE |
---|
1025 | CALL cv3_mixscale(nloc, ncum, nd, ment, m) |
---|
1026 | IF (debut) THEN |
---|
1027 | PRINT *, ' cv3_mixscale-> ' |
---|
1028 | END IF !(debut) THEN |
---|
1029 | END IF |
---|
1030 | END IF |
---|
1031 | |
---|
1032 | IF (iflag_con==4) THEN |
---|
1033 | if (prt_level >= 9) & |
---|
1034 | PRINT *, 'cva_driver -> cv_mixing' |
---|
1035 | CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis, & |
---|
1036 | ph, t, q, qs, u, v, h, lv, qnk, & |
---|
1037 | hp, tv, tvp, ep, clw, cbmf, & |
---|
1038 | m, ment, qent, uent, vent, nent, sigij, elij) |
---|
1039 | END IF |
---|
1040 | |
---|
1041 | IF (debut) THEN |
---|
1042 | PRINT *, ' cv_mixing ->' |
---|
1043 | END IF !(debut) THEN |
---|
1044 | ! do i = 1,nd |
---|
1045 | ! print*,'cv_mixing-> i,ment ',i,(ment(1,i,j),j=1,nd) |
---|
1046 | ! enddo |
---|
1047 | |
---|
1048 | ! ------------------------------------------------------------------- |
---|
1049 | ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS |
---|
1050 | ! ------------------------------------------------------------------- |
---|
1051 | IF (iflag_con==3) THEN |
---|
1052 | IF (debut) THEN |
---|
1053 | PRINT *, ' cva_driver -> cv3_unsat ' |
---|
1054 | END IF !(debut) THEN |
---|
1055 | |
---|
1056 | if (prt_level >= 9) & |
---|
1057 | PRINT *, 'cva_driver -> cv3_unsat' |
---|
1058 | CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, & ! na->nd |
---|
1059 | t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, & |
---|
1060 | th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, & |
---|
1061 | ep, sigp, clw, & |
---|
1062 | m, ment, elij, delt, plcl, coef_clos, & |
---|
1063 | mp, qp, up, vp, trap, wt, water, evap, fondue, ice, & |
---|
1064 | faci, b, sigd, & |
---|
1065 | wdtrainA, wdtrainM) ! RomP |
---|
1066 | ! |
---|
1067 | IF (prt_level >= 10) THEN |
---|
1068 | Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue ' |
---|
1069 | DO k = 1,nd |
---|
1070 | write (6, '(i4,5(1x,e13.6))'), & |
---|
1071 | k, mp(igout,k), water(igout,k), ice(igout,k), & |
---|
1072 | evap(igout,k), fondue(igout,k) |
---|
1073 | ENDDO |
---|
1074 | Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainM ' |
---|
1075 | DO k = 1,nd |
---|
1076 | write (6, '(i4,2(1x,e13.6))'), & |
---|
1077 | k, wdtrainA(igout,k), wdtrainM(igout,k) |
---|
1078 | ENDDO |
---|
1079 | ENDIF |
---|
1080 | ! |
---|
1081 | END IF !(iflag_con==3) |
---|
1082 | |
---|
1083 | IF (iflag_con==4) THEN |
---|
1084 | if (prt_level >= 9) & |
---|
1085 | PRINT *, 'cva_driver -> cv_unsat' |
---|
1086 | CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph, & |
---|
1087 | h, lv, ep, sigp, clw, m, ment, elij, & |
---|
1088 | iflag, mp, qp, up, vp, wt, water, evap) |
---|
1089 | END IF |
---|
1090 | |
---|
1091 | IF (debut) THEN |
---|
1092 | PRINT *, 'cv_unsat-> ' |
---|
1093 | END IF !(debut) THEN |
---|
1094 | |
---|
1095 | ! print *,'cv_unsat-> mp ',mp |
---|
1096 | ! print *,'cv_unsat-> water ',water |
---|
1097 | ! ------------------------------------------------------------------- |
---|
1098 | ! --- YIELD |
---|
1099 | ! (tendencies, precipitation, variables of interface with other |
---|
1100 | ! processes, etc) |
---|
1101 | ! ------------------------------------------------------------------- |
---|
1102 | |
---|
1103 | IF (iflag_con==3) THEN |
---|
1104 | |
---|
1105 | if (prt_level >= 9) & |
---|
1106 | PRINT *, 'cva_driver -> cv3_yield' |
---|
1107 | CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, & ! na->nd |
---|
1108 | icb, inb, delt, & |
---|
1109 | t, q, t_wake, q_wake, s_wake, u, v, tra, & |
---|
1110 | gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, & |
---|
1111 | ep, clw, m, tp, mp, qp, up, vp, trap, & |
---|
1112 | wt, water, ice, evap, fondue, faci, b, sigd, & |
---|
1113 | ment, qent, hent, iflag_mix, uent, vent, & |
---|
1114 | nent, elij, traent, sig, & |
---|
1115 | tv, tvp, wghti, & |
---|
1116 | iflag, precip, Vprecip, Vprecipi, ft, fq, fu, fv, ftra, & ! jyg |
---|
1117 | cbmf, upwd, dnwd, dnwd0, ma, mip, & |
---|
1118 | !! tls, tps, & ! useless . jyg |
---|
1119 | qcondc, wd, & |
---|
1120 | ftd, fqd, qnk, qtc, sigt, tau_cld_cv, coefw_cld_cv) |
---|
1121 | ! |
---|
1122 | IF (debut) THEN |
---|
1123 | PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1) |
---|
1124 | END IF !(debut) THEN |
---|
1125 | ! |
---|
1126 | IF (prt_level >= 10) THEN |
---|
1127 | Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', & |
---|
1128 | ft(igout,1), ftd(igout,1) |
---|
1129 | Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', & |
---|
1130 | fq(igout,1), fqd(igout,1) |
---|
1131 | ENDIF |
---|
1132 | ! |
---|
1133 | END IF |
---|
1134 | |
---|
1135 | IF (iflag_con==4) THEN |
---|
1136 | if (prt_level >= 9) & |
---|
1137 | PRINT *, 'cva_driver -> cv_yield' |
---|
1138 | CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt, & |
---|
1139 | t, q, u, v, & |
---|
1140 | gz, p, ph, h, hp, lv, cpn, & |
---|
1141 | ep, clw, frac, m, mp, qp, up, vp, & |
---|
1142 | wt, water, evap, & |
---|
1143 | ment, qent, uent, vent, nent, elij, & |
---|
1144 | tv, tvp, & |
---|
1145 | iflag, wd, qprime, tprime, & |
---|
1146 | precip, cbmf, ft, fq, fu, fv, ma, qcondc) |
---|
1147 | END IF |
---|
1148 | |
---|
1149 | !AC! |
---|
1150 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
1151 | !--- passive tracers |
---|
1152 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
1153 | |
---|
1154 | IF (iflag_con==3) THEN |
---|
1155 | !RomP >>> |
---|
1156 | if (prt_level >= 9) & |
---|
1157 | PRINT *, 'cva_driver -> cv3_tracer' |
---|
1158 | CALL cv3_tracer(nloc, len, ncum, nd, nd, & |
---|
1159 | ment, sigij, da, phi, phi2, d1a, dam, & |
---|
1160 | ep, vprecip, elij, clw, epmlmMm, eplaMm, & |
---|
1161 | icb, inb) |
---|
1162 | !RomP <<< |
---|
1163 | END IF |
---|
1164 | |
---|
1165 | !AC! |
---|
1166 | |
---|
1167 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
1168 | ! --- UNCOMPRESS THE FIELDS |
---|
1169 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
1170 | |
---|
1171 | |
---|
1172 | IF (iflag_con==3) THEN |
---|
1173 | if (prt_level >= 9) & |
---|
1174 | PRINT *, 'cva_driver -> cv3a_uncompress' |
---|
1175 | CALL cv3a_uncompress(nloc, len, ncum, nd, ntra, idcum, compress, & |
---|
1176 | iflag, icb, inb, & |
---|
1177 | precip, cbmf, plcl, plfc, wbeff, sig, w0, ptop2, & |
---|
1178 | ft, fq, fu, fv, ftra, & |
---|
1179 | sigd, ma, mip, vprecip, vprecipi, upwd, dnwd, dnwd0, & |
---|
1180 | qcondc, wd, cape, cin, & |
---|
1181 | tvp, & |
---|
1182 | ftd, fqd, & |
---|
1183 | Plim1, plim2, asupmax, supmax0, & |
---|
1184 | asupmaxmin, & |
---|
1185 | da, phi, mp, phi2, d1a, dam, sigij, & ! RomP |
---|
1186 | clw, elij, evap, ep, epmlmMm, eplaMm, & ! RomP |
---|
1187 | wdtrainA, wdtrainM, & ! RomP |
---|
1188 | qtc, sigt, epmax_diag, & ! epmax_cape |
---|
1189 | iflag1, kbas1, ktop1, & |
---|
1190 | precip1, cbmf1, plcl1, plfc1, wbeff1, sig1, w01, ptop21, & |
---|
1191 | ft1, fq1, fu1, fv1, ftra1, & |
---|
1192 | sigd1, ma1, mip1, vprecip1, vprecipi1, upwd1, dnwd1, dnwd01, & |
---|
1193 | qcondc1, wd1, cape1, cin1, & |
---|
1194 | tvp1, & |
---|
1195 | ftd1, fqd1, & |
---|
1196 | Plim11, plim21, asupmax1, supmax01, & |
---|
1197 | asupmaxmin1, & |
---|
1198 | da1, phi1, mp1, phi21, d1a1, dam1, sigij1, & ! RomP |
---|
1199 | clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, & ! RomP |
---|
1200 | wdtrainA1, wdtrainM1, & ! RomP |
---|
1201 | qtc1, sigt1, epmax_diag1) ! epmax_cape |
---|
1202 | END IF |
---|
1203 | |
---|
1204 | IF (iflag_con==4) THEN |
---|
1205 | if (prt_level >= 9) & |
---|
1206 | PRINT *, 'cva_driver -> cv_uncompress' |
---|
1207 | CALL cv_uncompress(nloc, len, ncum, nd, idcum, & |
---|
1208 | iflag, & |
---|
1209 | precip, cbmf, & |
---|
1210 | ft, fq, fu, fv, & |
---|
1211 | ma, qcondc, & |
---|
1212 | iflag1, & |
---|
1213 | precip1,cbmf1, & |
---|
1214 | ft1, fq1, fu1, fv1, & |
---|
1215 | ma1, qcondc1) |
---|
1216 | END IF |
---|
1217 | |
---|
1218 | END IF ! ncum>0 |
---|
1219 | |
---|
1220 | ! |
---|
1221 | ! In order take into account the possibility of changing the compression, |
---|
1222 | ! reset m, sig and w0 to zero for non-convective points. |
---|
1223 | DO k = 1,nd-1 |
---|
1224 | sig1(:, k) = sig1(:, k)*coef_convective(:) |
---|
1225 | w01(:, k) = w01(:, k)*coef_convective(:) |
---|
1226 | ENDDO |
---|
1227 | |
---|
1228 | IF (debut) THEN |
---|
1229 | PRINT *, ' cv_uncompress -> ' |
---|
1230 | debut = .FALSE. |
---|
1231 | END IF !(debut) THEN |
---|
1232 | |
---|
1233 | |
---|
1234 | RETURN |
---|
1235 | END SUBROUTINE cva_driver |
---|