1 | module cv3a_driver_mod |
---|
2 | contains |
---|
3 | |
---|
4 | ! ------------------------------------------------------------------- |
---|
5 | ! Prolog by Kerry Emanuel. |
---|
6 | ! ------------------------------------------------------------------- |
---|
7 | ! --- ARGUMENTS |
---|
8 | ! ------------------------------------------------------------------- |
---|
9 | ! --- On input: |
---|
10 | ! |
---|
11 | ! t: Array of absolute temperature (K) of dimension ND, with first |
---|
12 | ! index corresponding to lowest model level. Note that this array |
---|
13 | ! will be altered by the subroutine if dry convective adjustment |
---|
14 | ! occurs and if IPBL is not equal to 0. |
---|
15 | ! |
---|
16 | ! q: Array of specific humidity (gm/gm) of dimension ND, with first |
---|
17 | ! index corresponding to lowest model level. Must be defined |
---|
18 | ! at same grid levels as T. Note that this array will be altered |
---|
19 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
20 | ! |
---|
21 | ! qs: Array of saturation specific humidity of dimension ND, with first |
---|
22 | ! index corresponding to lowest model level. Must be defined |
---|
23 | ! at same grid levels as T. Note that this array will be altered |
---|
24 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
25 | ! |
---|
26 | ! t_wake: Array of absolute temperature (K), seen by unsaturated draughts, |
---|
27 | ! of dimension ND, with first index corresponding to lowest model level. |
---|
28 | ! |
---|
29 | ! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts, |
---|
30 | ! of dimension ND, with first index corresponding to lowest model level. |
---|
31 | ! Must be defined at same grid levels as T. |
---|
32 | ! |
---|
33 | ! qs_wake: Array of saturation specific humidity, seen by unsaturated draughts, |
---|
34 | ! of dimension ND, with first index corresponding to lowest model level. |
---|
35 | ! Must be defined at same grid levels as T. |
---|
36 | ! |
---|
37 | ! s_wake: Array of fractionnal area occupied by the wakes. |
---|
38 | ! |
---|
39 | ! u: Array of zonal wind velocity (m/s) of dimension ND, witth first |
---|
40 | ! index corresponding with the lowest model level. Defined at |
---|
41 | ! same levels as T. Note that this array will be altered if |
---|
42 | ! dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
43 | ! |
---|
44 | ! v: Same as u but for meridional velocity. |
---|
45 | ! |
---|
46 | ! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA), |
---|
47 | ! where NTRA is the number of different tracers. If no |
---|
48 | ! convective tracer transport is needed, define a dummy |
---|
49 | ! input array of dimension (ND,1). Tracers are defined at |
---|
50 | ! same vertical levels as T. Note that this array will be altered |
---|
51 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
52 | ! |
---|
53 | ! p: Array of pressure (mb) of dimension ND, with first |
---|
54 | ! index corresponding to lowest model level. Must be defined |
---|
55 | ! at same grid levels as T. |
---|
56 | ! |
---|
57 | ! ph: Array of pressure (mb) of dimension ND+1, with first index |
---|
58 | ! corresponding to lowest level. These pressures are defined at |
---|
59 | ! levels intermediate between those of P, T, Q and QS. The first |
---|
60 | ! value of PH should be greater than (i.e. at a lower level than) |
---|
61 | ! the first value of the array P. |
---|
62 | ! |
---|
63 | ! ALE: Available lifting Energy |
---|
64 | ! |
---|
65 | ! ALP: Available lifting Power |
---|
66 | ! |
---|
67 | ! nl: The maximum number of levels to which convection can penetrate, plus 1. |
---|
68 | ! NL MUST be less than or equal to ND-1. |
---|
69 | ! |
---|
70 | ! delt: The model time step (sec) between calls to CONVECT |
---|
71 | ! |
---|
72 | ! ---------------------------------------------------------------------------- |
---|
73 | ! --- On Output: |
---|
74 | ! |
---|
75 | ! iflag: An output integer whose value denotes the following: |
---|
76 | ! VALUE INTERPRETATION |
---|
77 | ! ----- -------------- |
---|
78 | ! 0 Moist convection occurs. |
---|
79 | ! 1 Moist convection occurs, but a CFL condition |
---|
80 | ! on the subsidence warming is violated. This |
---|
81 | ! does not cause the scheme to terminate. |
---|
82 | ! 2 Moist convection, but no precip because ep(inb) lt 0.0001 |
---|
83 | ! 3 No moist convection because new cbmf is 0 and old cbmf is 0. |
---|
84 | ! 4 No moist convection; atmosphere is not |
---|
85 | ! unstable |
---|
86 | ! 6 No moist convection because ihmin le minorig. |
---|
87 | ! 7 No moist convection because unreasonable |
---|
88 | ! parcel level temperature or specific humidity. |
---|
89 | ! 8 No moist convection: lifted condensation |
---|
90 | ! level is above the 200 mb level. |
---|
91 | ! 9 No moist convection: cloud base is higher |
---|
92 | ! then the level NL-1. |
---|
93 | ! 10 No moist convection: cloud top is too warm. |
---|
94 | ! |
---|
95 | ! |
---|
96 | ! ft: Array of temperature tendency (K/s) of dimension ND, defined at same |
---|
97 | ! grid levels as T, Q, QS and P. |
---|
98 | ! |
---|
99 | ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND, |
---|
100 | ! defined at same grid levels as T, Q, QS and P. |
---|
101 | ! |
---|
102 | ! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND, |
---|
103 | ! defined at same grid levels as T. |
---|
104 | ! |
---|
105 | ! fv: Same as FU, but for forcing of meridional velocity. |
---|
106 | ! |
---|
107 | ! ftra: Array of forcing of tracer content, in tracer mixing ratio per |
---|
108 | ! second, defined at same levels as T. Dimensioned (ND,NTRA). |
---|
109 | ! |
---|
110 | ! precip: Scalar convective precipitation rate (mm/day). |
---|
111 | ! |
---|
112 | ! wd: A convective downdraft velocity scale. For use in surface |
---|
113 | ! flux parameterizations. See convect.ps file for details. |
---|
114 | ! |
---|
115 | ! tprime: A convective downdraft temperature perturbation scale (K). |
---|
116 | ! For use in surface flux parameterizations. See convect.ps |
---|
117 | ! file for details. |
---|
118 | ! |
---|
119 | ! qprime: A convective downdraft specific humidity |
---|
120 | ! perturbation scale (gm/gm). |
---|
121 | ! For use in surface flux parameterizations. See convect.ps |
---|
122 | ! file for details. |
---|
123 | |
---|
124 | ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST |
---|
125 | ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT |
---|
126 | ! ITS NEXT CALL. That is, the value of CBMF must be "remembered" |
---|
127 | ! by the calling program between calls to CONVECT. |
---|
128 | ! |
---|
129 | ! det: Array of detrainment mass flux of dimension ND. |
---|
130 | ! ------------------------------------------------------------------- |
---|
131 | ! written by : Sandrine Bony-Lena , 17/05/2003, 11.19.41 |
---|
132 | ! S. Bony, Mar 2002: |
---|
133 | ! * Several modules corresponding to different physical processes |
---|
134 | ! * Several versions of convect may be used: |
---|
135 | ! - iflag_con=3: version lmd (previously named convect3) |
---|
136 | ! - iflag_con=4: version 4.3b (vect. version, previously convect1/2) |
---|
137 | ! + tard: - iflag_con=5: version lmd with ice (previously named convectg) |
---|
138 | ! S. Bony, Oct 2002: |
---|
139 | ! * Vectorization of convect3 (ie version lmd) |
---|
140 | SUBROUTINE cv3a_driver(len, nd, ndp1, ntra, nloc, k_upper, & |
---|
141 | iflag_con, iflag_mix, iflag_ice_thermo, iflag_clos, ok_conserv_q, & |
---|
142 | delt, & |
---|
143 | t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, & |
---|
144 | u1, v1, & |
---|
145 | p1, ph1, & |
---|
146 | Ale1, Alp1, omega1, & |
---|
147 | sig1feed1, sig2feed1, wght1, & |
---|
148 | iflag1, ft1, fq1, fu1, fv1, ftra1, & |
---|
149 | precip1, kbas1, ktop1, & |
---|
150 | cbmf1, plcl1, plfc1, wbeff1, & |
---|
151 | sig1, w01, & !input/output |
---|
152 | ptop21, sigd1, & |
---|
153 | ma1, mip1, Vprecip1, Vprecipi1, upwd1, dnwd1, dnwd01, & |
---|
154 | qcondc1, wd1, & |
---|
155 | cape1, cin1, tvp1, & |
---|
156 | ftd1, fqd1, & |
---|
157 | Plim11, Plim21, asupmax1, supmax01, asupmaxmin1, & |
---|
158 | da1, phi1, mp1, phi21, d1a1, dam1, sigij1, wghti1, & |
---|
159 | qta1, clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, & |
---|
160 | wdtrainA1, wdtrainS1, wdtrainM1, qtc1, sigt1, tau_cld_cv, & |
---|
161 | coefw_cld_cv, & |
---|
162 | epmax_diag1) |
---|
163 | |
---|
164 | USE print_control_mod, ONLY: prt_level, lunout |
---|
165 | USE cv3a_compress_mod |
---|
166 | USE cv3p_mixing_mod, ONLY : cv3p_mixing |
---|
167 | IMPLICIT NONE |
---|
168 | ! Input |
---|
169 | INTEGER, INTENT(IN) :: len ! first (i) dimension |
---|
170 | INTEGER, INTENT(IN) :: nd ! vertical (k) dimension |
---|
171 | INTEGER, INTENT(IN) :: ndp1 ! nd + 1 |
---|
172 | INTEGER, INTENT(IN) :: ntra ! number of tracors |
---|
173 | INTEGER, INTENT(IN) :: nloc ! dimension of arrays for compressed fields (nloc=len) pour l'instant |
---|
174 | INTEGER, INTENT(IN) :: k_upper ! upmost level for vertical loops |
---|
175 | INTEGER, INTENT(IN) :: iflag_con ! version of convect (3) |
---|
176 | INTEGER, INTENT(IN) :: iflag_mix ! version of mixing (0/1/2) |
---|
177 | INTEGER, INTENT(IN) :: iflag_ice_thermo ! accounting for ice thermodynamics (0/1) |
---|
178 | INTEGER, INTENT(IN) :: iflag_clos ! version of closure (0/1) |
---|
179 | LOGICAL, INTENT(IN) :: ok_conserv_q ! when true corrections for water conservation are swtiched on |
---|
180 | REAL, INTENT(IN) :: tau_cld_cv ! characteristic time of dissipation of mixing fluxes |
---|
181 | REAL, INTENT(IN) :: coefw_cld_cv ! coefficient for updraft velocity in convection |
---|
182 | REAL, INTENT(IN) :: delt ! time step |
---|
183 | REAL, DIMENSION(len, nd), INTENT(IN) :: t1 ! temperature (sat draught envt) |
---|
184 | REAL, DIMENSION(len, nd), INTENT(IN) :: q1 ! specific hum (sat draught envt) |
---|
185 | REAL, DIMENSION(len, nd), INTENT(IN) :: qs1 ! sat specific hum (sat draught envt) |
---|
186 | REAL, DIMENSION(len, nd), INTENT(IN) :: t1_wake ! temperature (unsat draught envt) |
---|
187 | REAL, DIMENSION(len, nd), INTENT(IN) :: q1_wake ! specific hum(unsat draught envt) |
---|
188 | REAL, DIMENSION(len, nd), INTENT(IN) :: qs1_wake ! sat specific hum(unsat draughts envt) |
---|
189 | REAL, DIMENSION(len), INTENT(IN) :: s1_wake ! fractionnal area covered by wakes |
---|
190 | REAL, DIMENSION(len, nd), INTENT(IN) :: u1 ! u-wind |
---|
191 | REAL, DIMENSION(len, nd), INTENT(IN) :: v1 ! v-wind |
---|
192 | REAL, DIMENSION(len, nd), INTENT(IN) :: p1 ! full level pressure |
---|
193 | REAL, DIMENSION(len, ndp1), INTENT(IN) :: ph1 ! half level pressure |
---|
194 | REAL, DIMENSION(len), INTENT(IN) :: Ale1 ! Available lifting Energy |
---|
195 | REAL, DIMENSION(len), INTENT(IN) :: Alp1 ! Available lifting Power |
---|
196 | REAL, DIMENSION(len, nd), INTENT(IN) :: omega1 |
---|
197 | REAL, INTENT(IN) :: sig1feed1 ! sigma coord/pressure at lower bound of feeding layer |
---|
198 | REAL, INTENT(IN) :: sig2feed1 ! sigma coord/pressure at upper bound of feeding layer |
---|
199 | REAL, DIMENSION(nd), INTENT(IN) :: wght1 ! weight density determining the feeding mixture |
---|
200 | |
---|
201 | ! Input/Output |
---|
202 | REAL, DIMENSION(len, nd), INTENT(INOUT) :: sig1 ! section adiabatic updraft |
---|
203 | REAL, DIMENSION(len, nd), INTENT(INOUT) :: w01 ! vertical velocity within adiab updraft |
---|
204 | |
---|
205 | ! Output |
---|
206 | INTEGER, DIMENSION(len), INTENT(OUT) :: iflag1 ! flag for Emanuel conditions |
---|
207 | REAL, DIMENSION(len, nd), INTENT(OUT) :: ft1 ! temp tend |
---|
208 | REAL, DIMENSION(len, nd), INTENT(OUT) :: fq1 ! spec hum tend |
---|
209 | REAL, DIMENSION(len, nd), INTENT(OUT) :: fu1 ! u-wind tend |
---|
210 | REAL, DIMENSION(len, nd), INTENT(OUT) :: fv1 ! v-wind tend |
---|
211 | REAL, DIMENSION(len, nd, ntra), INTENT(OUT) :: ftra1 ! tracor tend |
---|
212 | REAL, DIMENSION(len), INTENT(OUT) :: precip1 ! precipitation |
---|
213 | INTEGER, DIMENSION(len), INTENT(OUT) :: kbas1 ! cloud base level |
---|
214 | INTEGER, DIMENSION(len), INTENT(OUT) :: ktop1 ! cloud top level |
---|
215 | REAL, DIMENSION(len), INTENT(OUT) :: cbmf1 ! cloud base mass flux |
---|
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 ! top of entraining zone |
---|
220 | REAL, DIMENSION(len), INTENT(OUT) :: sigd1 |
---|
221 | REAL, DIMENSION(len, nd), INTENT(OUT) :: ma1 ! mass flux adiabatic updraft (staggered grid) |
---|
222 | REAL, DIMENSION(len, nd), INTENT(OUT) :: mip1 ! mass flux shed by the adiabatic updraft (extensive) |
---|
223 | REAL, DIMENSION(len, ndp1), INTENT(OUT) :: vprecip1 ! vertical profile of total precipitation (staggered grid) |
---|
224 | REAL, DIMENSION(len, ndp1), INTENT(OUT) :: vprecipi1 ! vertical profile of ice precipitation (staggered grid) |
---|
225 | REAL, DIMENSION(len, nd), INTENT(OUT) :: upwd1 ! total upward mass flux (adiab+mixed) (staggered grid) |
---|
226 | REAL, DIMENSION(len, nd), INTENT(OUT) :: dnwd1 ! saturated downward mass flux (mixed) (staggered grid) |
---|
227 | REAL, DIMENSION(len, nd), INTENT(OUT) :: dnwd01 ! unsaturated downward mass flux (staggered grid) |
---|
228 | REAL, DIMENSION(len, nd), INTENT(OUT) :: qcondc1 ! in-cld mixing ratio of condensed water (intensive) |
---|
229 | REAL, DIMENSION(len), INTENT(OUT) :: wd1 ! downdraft velocity scale for sfc fluxes |
---|
230 | REAL, DIMENSION(len), INTENT(OUT) :: cape1 |
---|
231 | REAL, DIMENSION(len), INTENT(OUT) :: cin1 |
---|
232 | REAL, DIMENSION(len, nd), INTENT(OUT) :: tvp1 ! adiab lifted parcell virt temp |
---|
233 | REAL, DIMENSION(len, nd), INTENT(OUT) :: ftd1 ! temperature tendency due to precipitations (K/s) |
---|
234 | REAL, DIMENSION(len, nd), INTENT(OUT) :: fqd1 ! specific humidity tendencies due to precipitations ((gm/gm)/s) |
---|
235 | REAL, DIMENSION(len), INTENT(OUT) :: Plim11 |
---|
236 | REAL, DIMENSION(len), INTENT(OUT) :: Plim21 |
---|
237 | REAL, DIMENSION(len, nd), INTENT(OUT) :: asupmax1 ! Highest mixing fraction of mixed updraughts |
---|
238 | REAL, DIMENSION(len), INTENT(OUT) :: supmax01 |
---|
239 | REAL, DIMENSION(len), INTENT(OUT) :: asupmaxmin1 |
---|
240 | REAL, DIMENSION(len, nd), INTENT(OUT) :: qtc1 ! in cloud water content / specific humidity in convection (intensive) |
---|
241 | REAL, DIMENSION(len, nd), INTENT(OUT) :: sigt1 ! surface fraction in adiabatic updrafts / fract. cloud area (intensive) |
---|
242 | REAL, DIMENSION(len, nd), INTENT(OUT) :: wdtrainA1 ! precipitation ejected from adiabatic draught |
---|
243 | REAL, DIMENSION(len, nd), INTENT(OUT) :: wdtrainS1 ! precipitation detrained from shedding of adiabatic draught |
---|
244 | REAL, DIMENSION(len, nd), INTENT(OUT) :: wdtrainM1 ! precipitation detrained from mixed draughts |
---|
245 | REAL, DIMENSION(len, nd), INTENT(OUT) :: mp1 ! unsat. mass flux (staggered grid) |
---|
246 | REAL, DIMENSION(len, nd), INTENT(OUT) :: da1 ! detrained mass flux of adiab. asc. air (extensive) |
---|
247 | REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: phi1 ! mass flux of envt. air in mixed draughts (extensive) |
---|
248 | REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: epmlmMm1 ! (extensive) |
---|
249 | REAL, DIMENSION(len, nd), INTENT(OUT) :: eplaMm1 ! (extensive) |
---|
250 | REAL, DIMENSION(len, nd), INTENT(OUT) :: evap1 ! evaporation rate in precip. downdraft. (intensive) |
---|
251 | REAL, DIMENSION(len, nd), INTENT(OUT) :: ep1 |
---|
252 | REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: sigij1 ! mass fraction of env. air in mixed draughts (intensive) |
---|
253 | REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: elij1! cond. water per unit mass of mixed draughts (intensive) |
---|
254 | REAL, DIMENSION(len, nd), INTENT(OUT) :: qta1 ! total water per unit mass of the adiab. asc. (intensive) |
---|
255 | REAL, DIMENSION(len, nd), INTENT(OUT) :: clw1 ! condensed water content of the adiabatic updraught / cond. water per unit mass of the adiab. asc. (intensive) |
---|
256 | REAL, DIMENSION(len, nd), INTENT(OUT) :: wghti1 ! final weight of the feeding layers (extensive) |
---|
257 | REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: phi21 ! (extensive) |
---|
258 | REAL, DIMENSION(len, nd), INTENT(OUT) :: d1a1 ! (extensive) |
---|
259 | REAL, DIMENSION(len, nd), INTENT(OUT) :: dam1 ! (extensive) |
---|
260 | REAL, DIMENSION(len), INTENT(OUT) :: epmax_diag1 |
---|
261 | |
---|
262 | ! Local (non compressed) arrays |
---|
263 | INTEGER i, k, il |
---|
264 | INTEGER icbmax |
---|
265 | INTEGER nk1(len) |
---|
266 | INTEGER icb1(len) |
---|
267 | INTEGER icbs1(len) |
---|
268 | |
---|
269 | LOGICAL ok_inhib ! True => possible inhibition of convection by dryness |
---|
270 | LOGICAL, SAVE :: debut = .TRUE. |
---|
271 | !$omp THREADPRIVATE(debut) |
---|
272 | |
---|
273 | REAL tnk1(len) |
---|
274 | REAL thnk1(len) |
---|
275 | REAL qnk1(len) |
---|
276 | REAL gznk1(len) |
---|
277 | REAL qsnk1(len) |
---|
278 | REAL unk1(len) |
---|
279 | REAL vnk1(len) |
---|
280 | REAL cpnk1(len) |
---|
281 | REAL hnk1(len) |
---|
282 | REAL pbase1(len) |
---|
283 | REAL buoybase1(len) |
---|
284 | |
---|
285 | REAL lf1(len, nd), lf1_wake(len, nd) |
---|
286 | REAL lv1(len, nd), lv1_wake(len, nd) |
---|
287 | REAL cpn1(len, nd), cpn1_wake(len, nd) |
---|
288 | REAL tv1(len, nd), tv1_wake(len, nd) |
---|
289 | REAL gz1(len, nd), gz1_wake(len, nd) |
---|
290 | REAL hm1(len, nd) |
---|
291 | REAL h1(len, nd), h1_wake(len, nd) |
---|
292 | REAL tp1(len, nd) |
---|
293 | REAL th1(len, nd), th1_wake(len, nd) |
---|
294 | |
---|
295 | REAL bid(len, nd) ! dummy array |
---|
296 | |
---|
297 | INTEGER ncum |
---|
298 | |
---|
299 | REAL p1feed1(len) ! pressure at lower bound of feeding layer |
---|
300 | REAL p2feed1(len) ! pressure at upper bound of feeding layer |
---|
301 | |
---|
302 | ! (local) compressed fields: |
---|
303 | LOGICAL compress ! True if compression occurs |
---|
304 | INTEGER iflag(nloc), nk(nloc), icb(nloc) |
---|
305 | INTEGER nent(nloc, nd) |
---|
306 | INTEGER icbs(nloc) |
---|
307 | INTEGER inb(nloc) |
---|
308 | |
---|
309 | REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc) |
---|
310 | REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd) |
---|
311 | REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd) |
---|
312 | REAL s_wake(nloc) |
---|
313 | REAL u(nloc, nd), v(nloc, nd) |
---|
314 | REAL gz(nloc, nd), h(nloc, nd) |
---|
315 | REAL h_wake(nloc, nd) |
---|
316 | REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd) |
---|
317 | REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd) |
---|
318 | REAL p(nloc, nd), ph(nloc, nd + 1), tv(nloc, nd), tp(nloc, nd) |
---|
319 | REAL tv_wake(nloc, nd) |
---|
320 | REAL clw(nloc, nd) |
---|
321 | REAL, DIMENSION(nloc, nd) :: qta, qpreca |
---|
322 | REAL pbase(nloc), buoybase(nloc), th(nloc, nd) |
---|
323 | REAL th_wake(nloc, nd) |
---|
324 | REAL tvp(nloc, nd) |
---|
325 | REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc) |
---|
326 | REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd) |
---|
327 | REAL buoy(nloc, nd) |
---|
328 | REAL cape(nloc) |
---|
329 | REAL cin(nloc) |
---|
330 | REAL m(nloc, nd) |
---|
331 | REAL mm(nloc, nd) |
---|
332 | REAL ment(nloc, nd, nd), sigij(nloc, nd, nd) |
---|
333 | REAL qent(nloc, nd, nd) |
---|
334 | REAL hent(nloc, nd, nd) |
---|
335 | REAL uent(nloc, nd, nd), vent(nloc, nd, nd) |
---|
336 | REAL ments(nloc, nd, nd), qents(nloc, nd, nd) |
---|
337 | REAL elij(nloc, nd, nd) |
---|
338 | REAL supmax(nloc, nd) |
---|
339 | REAL Ale(nloc), Alp(nloc), coef_clos(nloc) |
---|
340 | REAL omega(nloc, nd) |
---|
341 | REAL sigd(nloc) |
---|
342 | REAL, DIMENSION(len, nd) :: mp, qp, up, vp |
---|
343 | REAL, DIMENSION(len, nd) :: wt, water, evap |
---|
344 | REAL, DIMENSION(len, nd) :: ice, fondue, b |
---|
345 | REAL, DIMENSION(len, nd) :: frac_a, frac_s, faci |
---|
346 | REAL ft(nloc, nd), fq(nloc, nd) |
---|
347 | REAL ftd(nloc, nd), fqd(nloc, nd) |
---|
348 | REAL fu(nloc, nd), fv(nloc, nd) |
---|
349 | REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd) |
---|
350 | REAL ma(nloc, nd), mip(nloc, nd) |
---|
351 | REAL precip(nloc) |
---|
352 | REAL vprecip(nloc, nd + 1) |
---|
353 | REAL vprecipi(nloc, nd + 1) |
---|
354 | REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra) |
---|
355 | REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra) |
---|
356 | REAL qcondc(nloc, nd) |
---|
357 | REAL wd(nloc) |
---|
358 | REAL Plim1(nloc), plim2(nloc) |
---|
359 | REAL asupmax(nloc, nd) |
---|
360 | REAL supmax0(nloc) |
---|
361 | REAL asupmaxmin(nloc) |
---|
362 | |
---|
363 | REAL tnk(nloc), qnk(nloc), gznk(nloc) |
---|
364 | REAL wghti(nloc, nd) |
---|
365 | REAL hnk(nloc), unk(nloc), vnk(nloc) |
---|
366 | |
---|
367 | REAL qtc(nloc, nd) |
---|
368 | REAL sigt(nloc, nd) |
---|
369 | |
---|
370 | REAL wdtrainA(nloc, nd), wdtrainS(nloc, nd), wdtrainM(nloc, nd) |
---|
371 | REAL da(len, nd), phi(len, nd, nd) |
---|
372 | REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd) |
---|
373 | REAL phi2(len, nd, nd) |
---|
374 | REAL d1a(len, nd), dam(len, nd) |
---|
375 | REAL epmax_diag(nloc) ! epmax_cape |
---|
376 | |
---|
377 | CHARACTER(LEN=20) :: modname = 'cva_driver' |
---|
378 | CHARACTER(LEN=80) :: abort_message |
---|
379 | |
---|
380 | INTEGER, SAVE :: igout = 1 |
---|
381 | !$omp THREADPRIVATE(igout) |
---|
382 | |
---|
383 | type(compress_data_t) :: compress_data |
---|
384 | type(array_list) :: cv3a_compress_list, cv3a_uncompress_list |
---|
385 | |
---|
386 | if (iflag_con /= 3) call abort_physic("cv3a_driver", "iflag_con must be 3", 1) |
---|
387 | |
---|
388 | ! ------------------------------------------------------------------- |
---|
389 | ! --- SET CONSTANTS AND PARAMETERS |
---|
390 | ! ------------------------------------------------------------------- |
---|
391 | ! -- set simulation flags: |
---|
392 | ! (common cvflag |
---|
393 | CALL cv_flag(iflag_ice_thermo) |
---|
394 | ! -- set thermodynamical constants: |
---|
395 | ! (common cvthermo) |
---|
396 | CALL cv_thermo(iflag_con) |
---|
397 | ! -- set convect parameters |
---|
398 | ! includes microphysical parameters and parameters that |
---|
399 | ! control the rate of approach to quasi-equilibrium) |
---|
400 | ! (common cvparam) |
---|
401 | CALL cv3_param(nd, k_upper, delt) |
---|
402 | |
---|
403 | CALL cv3_incrcount(len, nd, delt, sig1) |
---|
404 | |
---|
405 | ! --------------------------------------------------------------------- |
---|
406 | ! --- INITIALIZE LOCAL ARRAYS AND PARAMETERS |
---|
407 | ! --------------------------------------------------------------------- |
---|
408 | |
---|
409 | DO il = 1, nloc |
---|
410 | coef_clos(il) = 1. |
---|
411 | END DO |
---|
412 | |
---|
413 | ! -------------------------------------------------------------------- |
---|
414 | ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY |
---|
415 | ! -------------------------------------------------------------------- |
---|
416 | |
---|
417 | IF (debut) PRINT *, 'Emanuel version 3 nouvelle' |
---|
418 | |
---|
419 | call driver_log('cv3_prelim') |
---|
420 | CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, & |
---|
421 | lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1) |
---|
422 | |
---|
423 | call driver_log('cv3_prelim') |
---|
424 | CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, & |
---|
425 | lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, & |
---|
426 | h1_wake, bid, th1_wake) |
---|
427 | |
---|
428 | ! -------------------------------------------------------------------- |
---|
429 | ! --- CONVECTIVE FEED |
---|
430 | ! -------------------------------------------------------------------- |
---|
431 | ! compute feeding layer potential temperature and mixing ratio : |
---|
432 | ! get bounds of feeding layer |
---|
433 | ! test niveaux couche alimentation KE |
---|
434 | IF (sig1feed1 == sig2feed1) THEN |
---|
435 | WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed' |
---|
436 | WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def' |
---|
437 | abort_message = '' |
---|
438 | CALL abort_physic(modname, abort_message, 1) |
---|
439 | END IF |
---|
440 | |
---|
441 | DO i = 1, len |
---|
442 | p1feed1(i) = sig1feed1*ph1(i, 1) |
---|
443 | p2feed1(i) = sig2feed1*ph1(i, 1) |
---|
444 | END DO |
---|
445 | |
---|
446 | call driver_log('cv3_feed') |
---|
447 | |
---|
448 | ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed |
---|
449 | iflag1(:) = 0 |
---|
450 | plcl1(:) = 0. |
---|
451 | wghti1(:, :) = 0. |
---|
452 | CALL cv3_feed(len, nd, ok_conserv_q, & |
---|
453 | t1, q1, u1, v1, p1, ph1, h1, gz1, & |
---|
454 | p1feed1, p2feed1, wght1, & |
---|
455 | wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, & |
---|
456 | cpnk1, hnk1, nk1, icb1, icbmax, iflag1, gznk1, plcl1) |
---|
457 | |
---|
458 | ! -------------------------------------------------------------------- |
---|
459 | ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part |
---|
460 | ! (up through ICB for convect4, up through ICB+1 for convect3) |
---|
461 | ! Calculates the lifted parcel virtual temperature at nk, the |
---|
462 | ! actual temperature, and the adiabatic liquid water content. |
---|
463 | ! -------------------------------------------------------------------- |
---|
464 | call driver_log('cv3_undilute1') |
---|
465 | ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed |
---|
466 | tvp1(:, :) = 0. |
---|
467 | clw1(:, :) = 0. |
---|
468 | CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, icb1, tnk1, qnk1, & ! nd->na |
---|
469 | gznk1, tp1, tvp1, clw1, icbs1) |
---|
470 | |
---|
471 | ! ------------------------------------------------------------------- |
---|
472 | ! --- TRIGGERING |
---|
473 | ! ------------------------------------------------------------------- |
---|
474 | |
---|
475 | call driver_log('cv3_trigger') |
---|
476 | CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na |
---|
477 | pbase1, buoybase1, iflag1, sig1, w01) |
---|
478 | |
---|
479 | ! ===================================================================== |
---|
480 | ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY |
---|
481 | ! ===================================================================== |
---|
482 | |
---|
483 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
484 | ! --- COMPRESS THE FIELDS |
---|
485 | ! (-> vectorization over convective gridpoints) |
---|
486 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
487 | compress = .true. |
---|
488 | if (compress) then |
---|
489 | compress_mode = COMPRESS_MODE_COMPRESS |
---|
490 | else |
---|
491 | compress_mode = COMPRESS_MODE_COPY |
---|
492 | endif |
---|
493 | |
---|
494 | call driver_log('cv3a_compress') |
---|
495 | |
---|
496 | call add_array_i1(cv3a_compress_list, iflag1, iflag) |
---|
497 | call add_array_i1(cv3a_compress_list, nk1, nk) |
---|
498 | call add_array_i1(cv3a_compress_list, icb1, icb) |
---|
499 | call add_array_i1(cv3a_compress_list, icbs1, icbs) |
---|
500 | call add_array_r1(cv3a_compress_list, plcl1, plcl) |
---|
501 | call add_array_r1(cv3a_compress_list, tnk1, tnk) |
---|
502 | call add_array_r1(cv3a_compress_list, qnk1, qnk) |
---|
503 | call add_array_r1(cv3a_compress_list, gznk1, gznk) |
---|
504 | call add_array_r1(cv3a_compress_list, hnk1, hnk) |
---|
505 | call add_array_r1(cv3a_compress_list, unk1, unk) |
---|
506 | call add_array_r1(cv3a_compress_list, vnk1, vnk) |
---|
507 | call add_array_r2(cv3a_compress_list, wghti1, wghti) |
---|
508 | call add_array_r1(cv3a_compress_list, pbase1, pbase) |
---|
509 | call add_array_r1(cv3a_compress_list, buoybase1, buoybase) |
---|
510 | call add_array_r2(cv3a_compress_list, th1, th) |
---|
511 | call add_array_r2(cv3a_compress_list, t1, t) |
---|
512 | call add_array_r2(cv3a_compress_list, q1, q) |
---|
513 | call add_array_r2(cv3a_compress_list, qs1, qs) |
---|
514 | call add_array_r2(cv3a_compress_list, t1_wake, t_wake) |
---|
515 | call add_array_r2(cv3a_compress_list, q1_wake, q_wake) |
---|
516 | call add_array_r2(cv3a_compress_list, qs1_wake, qs_wake) |
---|
517 | call add_array_r1(cv3a_compress_list, s1_wake, s_wake) |
---|
518 | call add_array_r2(cv3a_compress_list, u1, u) |
---|
519 | call add_array_r2(cv3a_compress_list, v1, v) |
---|
520 | call add_array_r2(cv3a_compress_list, gz1, gz) |
---|
521 | call add_array_r2(cv3a_compress_list, h1, h) |
---|
522 | call add_array_r2(cv3a_compress_list, th1_wake, th_wake) |
---|
523 | call add_array_r2(cv3a_compress_list, lv1, lv) |
---|
524 | call add_array_r2(cv3a_compress_list, lf1, lf) |
---|
525 | call add_array_r2(cv3a_compress_list, cpn1, cpn) |
---|
526 | call add_array_r2(cv3a_compress_list, p1, p) |
---|
527 | call add_array_r2(cv3a_compress_list, ph1, ph) |
---|
528 | call add_array_r2(cv3a_compress_list, tv1, tv) |
---|
529 | call add_array_r2(cv3a_compress_list, tp1, tp) |
---|
530 | call add_array_r2(cv3a_compress_list, tvp1, tvp) |
---|
531 | call add_array_r2(cv3a_compress_list, clw1, clw) |
---|
532 | call add_array_r2(cv3a_compress_list, h1_wake, h_wake) |
---|
533 | call add_array_r2(cv3a_compress_list, lv1_wake, lv_wake) |
---|
534 | call add_array_r2(cv3a_compress_list, lf1_wake, lf_wake) |
---|
535 | call add_array_r2(cv3a_compress_list, cpn1_wake, cpn_wake) |
---|
536 | call add_array_r2(cv3a_compress_list, tv1_wake, tv_wake) |
---|
537 | call add_array_r2(cv3a_compress_list, sig1, sig) |
---|
538 | call add_array_r1(cv3a_compress_list, sig1(:, nd), sig(:, nd)) |
---|
539 | call add_array_r2(cv3a_compress_list, w01, w0) |
---|
540 | call add_array_r1(cv3a_compress_list, Ale1, Ale) |
---|
541 | call add_array_r1(cv3a_compress_list, Alp1, Alp) |
---|
542 | call add_array_r2(cv3a_compress_list, omega1, omega) |
---|
543 | |
---|
544 | call cv3a_compress(len, (iflag1 == 0), cv3a_compress_list, compress_data) |
---|
545 | ncum = compress_data%ncum |
---|
546 | |
---|
547 | IF (ncum > 0) THEN |
---|
548 | |
---|
549 | ! ------------------------------------------------------------------- |
---|
550 | ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part : |
---|
551 | ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
---|
552 | ! --- & |
---|
553 | ! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE |
---|
554 | ! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD |
---|
555 | ! --- & |
---|
556 | ! --- FIND THE LEVEL OF NEUTRAL BUOYANCY |
---|
557 | ! ------------------------------------------------------------------- |
---|
558 | |
---|
559 | call driver_log('cv3_undilute2') |
---|
560 | CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, & |
---|
561 | tnk, qnk, gznk, hnk, t, q, qs, gz, & |
---|
562 | p, ph, h, tv, lv, lf, pbase, buoybase, plcl, & |
---|
563 | inb, tp, tvp, clw, hp, ep, sigp, buoy, & |
---|
564 | frac_a, frac_s, qpreca, qta) |
---|
565 | |
---|
566 | ! epmax_cape |
---|
567 | ! on recalcule ep et hp |
---|
568 | call driver_log('cv3_epmax_cape') |
---|
569 | call cv3_epmax_fn_cape(nloc, ncum, nd & |
---|
570 | , ep, hp, icb, inb, clw, nk, t, h, hnk, lv, lf, frac_s & |
---|
571 | , pbase, p, ph, tv, buoy, sig, w0, iflag & |
---|
572 | , epmax_diag) |
---|
573 | |
---|
574 | ! ------------------------------------------------------------------- |
---|
575 | ! --- MIXING(1) (if iflag_mix .ge. 1) |
---|
576 | ! ------------------------------------------------------------------- |
---|
577 | |
---|
578 | IF (iflag_mix >= 1) THEN |
---|
579 | CALL zilch(supmax, nloc*nd) |
---|
580 | call driver_log('cv3p_mixing') |
---|
581 | CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, & |
---|
582 | ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qta, & |
---|
583 | unk, vnk, hp, tv, tvp, ep, clw, sig, & |
---|
584 | ment, qent, hent, uent, vent, nent, & |
---|
585 | sigij, elij, supmax, ments, qents, traent) |
---|
586 | ELSE |
---|
587 | CALL zilch(supmax, nloc*nd) |
---|
588 | END IF |
---|
589 | |
---|
590 | ! ------------------------------------------------------------------- |
---|
591 | ! --- CLOSURE |
---|
592 | ! ------------------------------------------------------------------- |
---|
593 | |
---|
594 | ptop2(:) = 0 |
---|
595 | IF (iflag_clos == 0) THEN |
---|
596 | call driver_log('cv3_closure') |
---|
597 | CALL cv3_closure(nloc, ncum, nd, icb, inb, & |
---|
598 | pbase, p, ph, tv, buoy, & |
---|
599 | sig, w0, cape, m, iflag) |
---|
600 | END IF ! iflag_clos==0 |
---|
601 | |
---|
602 | ok_inhib = iflag_mix == 2 |
---|
603 | |
---|
604 | IF (iflag_clos == 1) PRINT *, ' pas d appel cv3p_closure' |
---|
605 | |
---|
606 | IF (iflag_clos == 2) THEN |
---|
607 | call driver_log('cv3p1_closure') |
---|
608 | CALL cv3p1_closure(nloc, ncum, nd, icb, inb, & |
---|
609 | pbase, plcl, p, ph, tv, tvp, buoy, & |
---|
610 | supmax, ok_inhib, Ale, Alp, omega, & |
---|
611 | sig, w0, ptop2, cape, cin, m, iflag, coef_clos, & |
---|
612 | Plim1, plim2, asupmax, supmax0, & |
---|
613 | asupmaxmin, cbmf, plfc, wbeff) |
---|
614 | if (prt_level >= 10) PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1) |
---|
615 | END IF ! iflag_clos==2 |
---|
616 | |
---|
617 | IF (iflag_clos == 3) THEN |
---|
618 | call driver_log('cv3p2_closure') |
---|
619 | CALL cv3p2_closure(nloc, ncum, nd, icb, inb, & |
---|
620 | pbase, plcl, p, ph, tv, tvp, buoy, & |
---|
621 | supmax, ok_inhib, Ale, Alp, omega, & |
---|
622 | sig, w0, ptop2, cape, cin, m, iflag, coef_clos, & |
---|
623 | Plim1, plim2, asupmax, supmax0, & |
---|
624 | asupmaxmin, cbmf, plfc, wbeff) |
---|
625 | if (prt_level >= 10) PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1) |
---|
626 | END IF ! iflag_clos==3 |
---|
627 | |
---|
628 | ! ------------------------------------------------------------------- |
---|
629 | ! --- MIXING(2) |
---|
630 | ! ------------------------------------------------------------------- |
---|
631 | |
---|
632 | IF (iflag_mix == 0) THEN |
---|
633 | call driver_log('cv3_mixing') |
---|
634 | CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, & ! na->nd |
---|
635 | ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qnk, & |
---|
636 | unk, vnk, hp, tv, tvp, ep, clw, m, sig, & |
---|
637 | ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent) |
---|
638 | CALL zilch(hent, nloc*nd*nd) |
---|
639 | ELSE |
---|
640 | mm(:, :) = m(:, :) |
---|
641 | CALL cv3_mixscale(nloc, ncum, nd, ment, mm) |
---|
642 | IF (debut) PRINT *, ' cv3_mixscale-> ' |
---|
643 | END IF |
---|
644 | |
---|
645 | IF (debut) PRINT *, ' cv_mixing ->' |
---|
646 | |
---|
647 | ! ------------------------------------------------------------------- |
---|
648 | ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS |
---|
649 | ! ------------------------------------------------------------------- |
---|
650 | IF (debut) PRINT *, ' cva_driver -> cv3_unsat ' |
---|
651 | |
---|
652 | call driver_log('cv3_unsat') |
---|
653 | CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, & |
---|
654 | t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, & |
---|
655 | th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, & |
---|
656 | ep, sigp, clw, frac_s, qpreca, frac_a, qta, & |
---|
657 | m, ment, elij, delt, plcl, coef_clos, & |
---|
658 | mp, qp, up, vp, trap, wt, water, evap, fondue, ice, & |
---|
659 | faci, b, sigd, & |
---|
660 | wdtrainA, wdtrainS, wdtrainM) |
---|
661 | IF (prt_level >= 10) THEN |
---|
662 | Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue ' |
---|
663 | DO k = 1, nd |
---|
664 | write (6, '(i4,5(1x,e13.6))'), & |
---|
665 | k, mp(igout, k), water(igout, k), ice(igout, k), & |
---|
666 | evap(igout, k), fondue(igout, k) |
---|
667 | ENDDO |
---|
668 | Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainS, wdtrainM ' |
---|
669 | DO k = 1, nd |
---|
670 | write (6, '(i4,3(1x,e13.6))'), & |
---|
671 | k, wdtrainA(igout, k), wdtrainS(igout, k), wdtrainM(igout, k) |
---|
672 | ENDDO |
---|
673 | ENDIF |
---|
674 | |
---|
675 | IF (debut) PRINT *, 'cv_unsat-> ' |
---|
676 | ! ------------------------------------------------------------------- |
---|
677 | ! YIELD |
---|
678 | ! (tendencies, precipitation, variables of interface with other processes, etc) |
---|
679 | ! ------------------------------------------------------------------- |
---|
680 | |
---|
681 | call driver_log('cv3_yield') |
---|
682 | CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, & |
---|
683 | icb, inb, delt, & |
---|
684 | t, q, t_wake, q_wake, s_wake, u, v, tra, & |
---|
685 | gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, & |
---|
686 | ep, clw, qpreca, m, tp, mp, qp, up, vp, trap, & |
---|
687 | wt, water, ice, evap, fondue, faci, b, sigd, & |
---|
688 | ment, qent, hent, iflag_mix, uent, vent, & |
---|
689 | nent, elij, traent, sig, & |
---|
690 | tv, tvp, wghti, & |
---|
691 | iflag, precip, Vprecip, Vprecipi, ft, fq, fu, fv, ftra, & |
---|
692 | cbmf, upwd, dnwd, dnwd0, ma, mip, & |
---|
693 | qcondc, wd, & |
---|
694 | ftd, fqd, qta, qtc, sigt, tau_cld_cv, coefw_cld_cv) |
---|
695 | |
---|
696 | ! Test conseravtion de l'eau |
---|
697 | IF (debut) PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1) |
---|
698 | IF (prt_level >= 10) THEN |
---|
699 | Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', & |
---|
700 | ft(igout, 1), ftd(igout, 1) |
---|
701 | Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', & |
---|
702 | fq(igout, 1), fqd(igout, 1) |
---|
703 | ENDIF |
---|
704 | |
---|
705 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
706 | !--- passive tracers |
---|
707 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
708 | |
---|
709 | call driver_log('cv3_tracer') |
---|
710 | CALL cv3_tracer(nloc, len, ncum, nd, nd, & |
---|
711 | ment, sigij, da, phi, phi2, d1a, dam, & |
---|
712 | ep, vprecip, elij, clw, epmlmMm, eplaMm, & |
---|
713 | icb, inb) |
---|
714 | END IF ! ncum>0 |
---|
715 | |
---|
716 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
717 | ! --- UNCOMPRESS THE FIELDS |
---|
718 | ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
719 | call driver_log('cv3a_uncompress') |
---|
720 | call add_array_i1(cv3a_uncompress_list, iflag, iflag1, init=.false.) |
---|
721 | call add_array_i1(cv3a_uncompress_list, icb, kbas1) |
---|
722 | call add_array_i1(cv3a_uncompress_list, inb, ktop1) |
---|
723 | call add_array_r1(cv3a_uncompress_list, precip, precip1) |
---|
724 | call add_array_r1(cv3a_uncompress_list, cbmf, cbmf1) |
---|
725 | call add_array_r1(cv3a_uncompress_list, plcl, plcl1, init=.false.) |
---|
726 | call add_array_r1(cv3a_uncompress_list, plfc, plfc1) |
---|
727 | call add_array_r1(cv3a_uncompress_list, wbeff, wbeff1) |
---|
728 | call add_array_r2(cv3a_uncompress_list, sig, sig1) |
---|
729 | call add_array_r2(cv3a_uncompress_list, w0, w01) |
---|
730 | call add_array_r1(cv3a_uncompress_list, ptop2, ptop21) |
---|
731 | call add_array_r2(cv3a_uncompress_list, ft, ft1) |
---|
732 | call add_array_r2(cv3a_uncompress_list, fq, fq1) |
---|
733 | call add_array_r2(cv3a_uncompress_list, fu, fu1) |
---|
734 | call add_array_r2(cv3a_uncompress_list, fv, fv1) |
---|
735 | call add_array_r1(cv3a_uncompress_list, sigd, sigd1) |
---|
736 | call add_array_r2(cv3a_uncompress_list, ma, ma1) |
---|
737 | call add_array_r2(cv3a_uncompress_list, mip, mip1) |
---|
738 | call add_array_r2(cv3a_uncompress_list, vprecip, vprecip1) |
---|
739 | call add_array_r2(cv3a_uncompress_list, vprecipi, vprecipi1) |
---|
740 | call add_array_r2(cv3a_uncompress_list, upwd, upwd1) |
---|
741 | call add_array_r2(cv3a_uncompress_list, dnwd, dnwd1) |
---|
742 | call add_array_r2(cv3a_uncompress_list, dnwd0, dnwd01) |
---|
743 | call add_array_r2(cv3a_uncompress_list, qcondc, qcondc1) |
---|
744 | call add_array_r1(cv3a_uncompress_list, wd, wd1) |
---|
745 | cape1(:) = -1. |
---|
746 | call add_array_r1(cv3a_uncompress_list, cape, cape1, init=.false.) |
---|
747 | cin1(:) = -100000. |
---|
748 | call add_array_r1(cv3a_uncompress_list, cin, cin1, init=.false.) |
---|
749 | call add_array_r2(cv3a_uncompress_list, tvp, tvp1, init=.false.) |
---|
750 | call add_array_r2(cv3a_uncompress_list, ftd, ftd1) |
---|
751 | call add_array_r2(cv3a_uncompress_list, fqd, fqd1) |
---|
752 | call add_array_r1(cv3a_uncompress_list, Plim1, Plim11) |
---|
753 | call add_array_r1(cv3a_uncompress_list, plim2, plim21) |
---|
754 | call add_array_r2(cv3a_uncompress_list, asupmax, asupmax1) |
---|
755 | call add_array_r1(cv3a_uncompress_list, supmax0, supmax01) |
---|
756 | call add_array_r1(cv3a_uncompress_list, asupmaxmin, asupmaxmin1) |
---|
757 | call add_array_r2(cv3a_uncompress_list, da, da1) |
---|
758 | call add_array_r3(cv3a_uncompress_list, phi, phi1) |
---|
759 | call add_array_r2(cv3a_uncompress_list, mp, mp1) |
---|
760 | call add_array_r3(cv3a_uncompress_list, phi2, phi21) |
---|
761 | call add_array_r2(cv3a_uncompress_list, d1a, d1a1) |
---|
762 | call add_array_r2(cv3a_uncompress_list, dam, dam1) |
---|
763 | call add_array_r3(cv3a_uncompress_list, sigij, sigij1) |
---|
764 | call add_array_r2(cv3a_uncompress_list, qta, qta1) |
---|
765 | call add_array_r2(cv3a_uncompress_list, clw, clw1, init=.false.) |
---|
766 | call add_array_r3(cv3a_uncompress_list, elij, elij1) |
---|
767 | call add_array_r2(cv3a_uncompress_list, evap, evap1) |
---|
768 | call add_array_r2(cv3a_uncompress_list, ep, ep1) |
---|
769 | call add_array_r3(cv3a_uncompress_list, epmlmMm, epmlmMm1) |
---|
770 | call add_array_r2(cv3a_uncompress_list, eplaMm, eplaMm1) |
---|
771 | call add_array_r2(cv3a_uncompress_list, wdtrainA, wdtrainA1) |
---|
772 | call add_array_r2(cv3a_uncompress_list, wdtrainS, wdtrainS1) |
---|
773 | call add_array_r2(cv3a_uncompress_list, wdtrainM, wdtrainM1) |
---|
774 | call add_array_r2(cv3a_uncompress_list, qtc, qtc1) |
---|
775 | call add_array_r2(cv3a_uncompress_list, sigt, sigt1) |
---|
776 | call add_array_r1(cv3a_uncompress_list, epmax_diag, epmax_diag1) |
---|
777 | call add_array_r1(cv3a_uncompress_list, sig(:, nd), sig1(:, nd)) |
---|
778 | call cv3a_uncompress(len, compress_data, cv3a_uncompress_list) |
---|
779 | |
---|
780 | IF (prt_level >= 10) THEN |
---|
781 | Print *, 'cva_driver after cv3_uncompress:ft1(1) , ftd1(1) ', & |
---|
782 | ft1(igout, 1), ftd1(igout, 1) |
---|
783 | Print *, 'cva_driver after cv3_uncompress:fq1(1) , fqd1(1) ', & |
---|
784 | fq1(igout, 1), fqd1(igout, 1) |
---|
785 | ENDIF |
---|
786 | IF (debut) THEN |
---|
787 | PRINT *, ' cv_uncompress -> ' |
---|
788 | debut = .FALSE. |
---|
789 | END IF |
---|
790 | |
---|
791 | ftra1(:,:,:) = 0 |
---|
792 | |
---|
793 | END SUBROUTINE |
---|
794 | |
---|
795 | subroutine driver_log(message) |
---|
796 | use print_control_mod, only: prt_level |
---|
797 | character(*) :: message |
---|
798 | if (prt_level >= 9) PRINT *, 'cva_driver ->', message |
---|
799 | end subroutine |
---|
800 | end module |
---|