/[lmdze]/trunk/phylmd/cv_driver.f
ViewVC logotype

Contents of /trunk/phylmd/cv_driver.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 99 - (show annotations)
Wed Jul 2 18:39:15 2014 UTC (9 years, 10 months ago) by guez
File size: 21486 byte(s)
Created procedure test_disvert (following LMDZ). Added procedures
hybrid and funcd in module disvert_m. Upgraded compute_ab from
internal procedure of disvert to module procedure. Added variables y,
ya in module disvert_m. Upgraded s from local variable of procedure
disvert to module variable.

Renamed allowed value of variable vert_sampling in procedure disvert
from "read" to "read_hybrid". Added possibility to read pressure
values, value "read_pressure". Replaced vertical distribution for
value "param" by the distribution "strato_correct" from LMDZ (but kept
the value "param"). In case "tropo", replaced 1 by dsigmin (following
LMDZ). In case "strato", replaced 0.3 by dsigmin (following LMDZ).

Changed computation of bp in procedure compute_ab.

Removed debugindex case in clmain. Removed useless argument rlon of
procedure clmain. Removed useless variables ytaux, ytauy of procedure
clmain.

Removed intermediary variables tsol, qsol, tsolsrf, tslab in procedure
etat0.

Removed variable ok_veget:. coupling with the model Orchid is not
possible. Removed variable ocean: modeling an ocean slab is not
possible.

Removed useless variables tmp_rriv and tmp_rcoa from module
interface_surf.

Moved initialization of variables da, mp, phi in procedure physiq to
to inside the test iflag_con >= 3.

1 module cv_driver_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv_driver(len, nd, t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, &
8 fq1, fu1, fv1, precip1, VPrecip1, cbmf1, sig1, w01, icb1, inb1, delt, &
9 Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
10
11 ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12 ! Main driver for convection
13 ! Author: S. Bony, March 2002
14
15 ! Several modules corresponding to different physical processes
16
17 ! Several versions of convect may be used:
18 ! - iflag_con = 3: version lmd
19 ! - iflag_con = 4: version 4.3b
20
21 use clesphys2, only: iflag_con
22 use cv3_compress_m, only: cv3_compress
23 use cv3_mixing_m, only: cv3_mixing
24 use cv3_param_m, only: cv3_param
25 use cv3_prelim_m, only: cv3_prelim
26 use cv3_tracer_m, only: cv3_tracer
27 use cv3_uncompress_m, only: cv3_uncompress
28 use cv3_unsat_m, only: cv3_unsat
29 use cv3_yield_m, only: cv3_yield
30 use cv_uncompress_m, only: cv_uncompress
31 USE dimphy, ONLY: klev, klon
32
33 integer, intent(in):: len ! first dimension
34 integer, intent(in):: nd ! vertical dimension
35 real, intent(in):: t1(len, nd) ! temperature
36 real q1(len, nd) ! Input specific hum
37 real qs1(len, nd)
38 ! qs1 Real Input sat specific hum
39 real, intent(in):: u1(len, nd)
40 ! u1 Real Input u-wind
41 real, intent(in):: v1(len, nd)
42 ! v1 Real Input v-wind
43 real p1(len, nd)
44 ! p1 Real Input full level pressure
45 real ph1(len, nd + 1)
46 ! ph1 Real Input half level pressure
47 integer iflag1(len)
48 ! iflag1 Integer Output flag for Emanuel conditions
49 real ft1(len, nd)
50 ! ft1 Real Output temp tend
51 real fq1(len, nd)
52 ! fq1 Real Output spec hum tend
53 real fu1(len, nd)
54 ! fu1 Real Output u-wind tend
55 real fv1(len, nd)
56 ! fv1 Real Output v-wind tend
57 real precip1(len)
58 ! precip1 Real Output precipitation
59 real VPrecip1(len, nd+1)
60 ! VPrecip1 Real Output vertical profile of precipitations
61 real cbmf1(len)
62 ! cbmf1 Real Output cloud base mass flux
63 real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft
64
65 real, intent(inout):: w01(klon, klev)
66 ! vertical velocity within adiabatic updraft
67
68 integer icb1(klon)
69 integer inb1(klon)
70 real, intent(in):: delt
71 ! delt Real Input time step
72 real Ma1(len, nd)
73 ! Ma1 Real Output mass flux adiabatic updraft
74 real, intent(out):: upwd1(len, nd) ! total upward mass flux (adiab+mixed)
75 real, intent(out):: dnwd1(len, nd) ! saturated downward mass flux (mixed)
76 real, intent(out):: dnwd01(len, nd) ! unsaturated downward mass flux
77
78 real qcondc1(len, nd) ! cld
79 ! qcondc1 Real Output in-cld mixing ratio of condensed water
80 real wd1(len) ! gust
81 ! wd1 Real Output downdraft velocity scale for sfc fluxes
82 real cape1(len)
83 ! cape1 Real Output CAPE
84
85 real, intent(inout):: da1(len, nd), phi1(len, nd, nd), mp1(len, nd)
86
87 !-------------------------------------------------------------------
88 ! --- ARGUMENTS
89 !-------------------------------------------------------------------
90 ! --- On input:
91
92 ! t: Array of absolute temperature (K) of dimension ND, with first
93 ! index corresponding to lowest model level. Note that this array
94 ! will be altered by the subroutine if dry convective adjustment
95 ! occurs and if IPBL is not equal to 0.
96
97 ! q: Array of specific humidity (gm/gm) of dimension ND, with first
98 ! index corresponding to lowest model level. Must be defined
99 ! at same grid levels as T. Note that this array will be altered
100 ! if dry convective adjustment occurs and if IPBL is not equal to 0.
101
102 ! qs: Array of saturation specific humidity of dimension ND, with first
103 ! index corresponding to lowest model level. Must be defined
104 ! at same grid levels as T. Note that this array will be altered
105 ! if dry convective adjustment occurs and if IPBL is not equal to 0.
106
107 ! u: Array of zonal wind velocity (m/s) of dimension ND, witth first
108 ! index corresponding with the lowest model level. Defined at
109 ! same levels as T. Note that this array will be altered if
110 ! dry convective adjustment occurs and if IPBL is not equal to 0.
111
112 ! v: Same as u but for meridional velocity.
113
114 ! p: Array of pressure (mb) of dimension ND, with first
115 ! index corresponding to lowest model level. Must be defined
116 ! at same grid levels as T.
117
118 ! ph: Array of pressure (mb) of dimension ND+1, with first index
119 ! corresponding to lowest level. These pressures are defined at
120 ! levels intermediate between those of P, T, Q and QS. The first
121 ! value of PH should be greater than (i.e. at a lower level than)
122 ! the first value of the array P.
123
124 ! nl: The maximum number of levels to which convection can penetrate, plus 1.
125 ! NL MUST be less than or equal to ND-1.
126
127 ! delt: The model time step (sec) between calls to CONVECT
128
129 !----------------------------------------------------------------------------
130 ! --- On Output:
131
132 ! iflag: An output integer whose value denotes the following:
133 ! VALUE INTERPRETATION
134 ! ----- --------------
135 ! 0 Moist convection occurs.
136 ! 1 Moist convection occurs, but a CFL condition
137 ! on the subsidence warming is violated. This
138 ! does not cause the scheme to terminate.
139 ! 2 Moist convection, but no precip because ep(inb) lt 0.0001
140 ! 3 No moist convection because new cbmf is 0 and old cbmf is 0.
141 ! 4 No moist convection; atmosphere is not
142 ! unstable
143 ! 6 No moist convection because ihmin le minorig.
144 ! 7 No moist convection because unreasonable
145 ! parcel level temperature or specific humidity.
146 ! 8 No moist convection: lifted condensation
147 ! level is above the 200 mb level.
148 ! 9 No moist convection: cloud base is higher
149 ! then the level NL-1.
150
151 ! ft: Array of temperature tendency (K/s) of dimension ND, defined at same
152 ! grid levels as T, Q, QS and P.
153
154 ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
155 ! defined at same grid levels as T, Q, QS and P.
156
157 ! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND,
158 ! defined at same grid levels as T.
159
160 ! fv: Same as FU, but for forcing of meridional velocity.
161
162 ! precip: Scalar convective precipitation rate (mm/day).
163
164 ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).
165
166 ! wd: A convective downdraft velocity scale. For use in surface
167 ! flux parameterizations. See convect.ps file for details.
168
169 ! tprime: A convective downdraft temperature perturbation scale (K).
170 ! For use in surface flux parameterizations. See convect.ps
171 ! file for details.
172
173 ! qprime: A convective downdraft specific humidity
174 ! perturbation scale (gm/gm).
175 ! For use in surface flux parameterizations. See convect.ps
176 ! file for details.
177
178 ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
179 ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
180 ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"
181 ! by the calling program between calls to CONVECT.
182
183 ! det: Array of detrainment mass flux of dimension ND.
184
185 !-------------------------------------------------------------------
186
187 ! Local arrays
188
189 real da(len, nd), phi(len, nd, nd), mp(len, nd)
190
191 integer i, k, il
192 integer icbmax
193 integer nk1(klon)
194 integer icbs1(klon)
195
196 real plcl1(klon)
197 real tnk1(klon)
198 real qnk1(klon)
199 real gznk1(klon)
200 real pbase1(klon)
201 real buoybase1(klon)
202
203 real lv1(klon, klev)
204 real cpn1(klon, klev)
205 real tv1(klon, klev)
206 real gz1(klon, klev)
207 real hm1(klon, klev)
208 real h1(klon, klev)
209 real tp1(klon, klev)
210 real tvp1(klon, klev)
211 real clw1(klon, klev)
212 real th1(klon, klev)
213
214 integer ncum
215
216 ! (local) compressed fields:
217
218 integer nloc
219 parameter (nloc = klon) ! pour l'instant
220
221 integer idcum(nloc)
222 integer iflag(nloc), nk(nloc), icb(nloc)
223 integer nent(nloc, klev)
224 integer icbs(nloc)
225 integer inb(nloc), inbis(nloc)
226
227 real cbmf(nloc), plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)
228 real t(nloc, klev), q(nloc, klev), qs(nloc, klev)
229 real u(nloc, klev), v(nloc, klev)
230 real gz(nloc, klev), h(nloc, klev), lv(nloc, klev), cpn(nloc, klev)
231 real p(nloc, klev), ph(nloc, klev+1), tv(nloc, klev), tp(nloc, klev)
232 real clw(nloc, klev)
233 real dph(nloc, klev)
234 real pbase(nloc), buoybase(nloc), th(nloc, klev)
235 real tvp(nloc, klev)
236 real sig(nloc, klev), w0(nloc, klev)
237 real hp(nloc, klev), ep(nloc, klev), sigp(nloc, klev)
238 real frac(nloc), buoy(nloc, klev)
239 real cape(nloc)
240 real m(nloc, klev), ment(nloc, klev, klev), qent(nloc, klev, klev)
241 real uent(nloc, klev, klev), vent(nloc, klev, klev)
242 real ments(nloc, klev, klev), qents(nloc, klev, klev)
243 real sij(nloc, klev, klev), elij(nloc, klev, klev)
244 real qp(nloc, klev), up(nloc, klev), vp(nloc, klev)
245 real wt(nloc, klev), water(nloc, klev), evap(nloc, klev)
246 real b(nloc, klev), ft(nloc, klev), fq(nloc, klev)
247 real fu(nloc, klev), fv(nloc, klev)
248 real upwd(nloc, klev), dnwd(nloc, klev), dnwd0(nloc, klev)
249 real Ma(nloc, klev), mike(nloc, klev), tls(nloc, klev)
250 real tps(nloc, klev), qprime(nloc), tprime(nloc)
251 real precip(nloc)
252 real VPrecip(nloc, klev+1)
253 real qcondc(nloc, klev) ! cld
254 real wd(nloc) ! gust
255
256 !-------------------------------------------------------------------
257 ! --- SET CONSTANTS AND PARAMETERS
258 !-------------------------------------------------------------------
259
260 ! -- set simulation flags:
261 ! (common cvflag)
262
263 CALL cv_flag
264
265 ! -- set thermodynamical constants:
266 ! (common cvthermo)
267
268 CALL cv_thermo
269
270 ! -- set convect parameters
271
272 ! includes microphysical parameters and parameters that
273 ! control the rate of approach to quasi-equilibrium)
274 ! (common cvparam)
275
276 if (iflag_con.eq.3) then
277 CALL cv3_param(nd, delt)
278 endif
279
280 if (iflag_con.eq.4) then
281 CALL cv_param(nd)
282 endif
283
284 !---------------------------------------------------------------------
285 ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
286 !---------------------------------------------------------------------
287
288 do k = 1, nd
289 do i = 1, len
290 ft1(i, k) = 0.0
291 fq1(i, k) = 0.0
292 fu1(i, k) = 0.0
293 fv1(i, k) = 0.0
294 tvp1(i, k) = 0.0
295 tp1(i, k) = 0.0
296 clw1(i, k) = 0.0
297 !ym
298 clw(i, k) = 0.0
299 gz1(i, k) = 0.
300 VPrecip1(i, k) = 0.
301 Ma1(i, k) = 0.0
302 upwd1(i, k) = 0.0
303 dnwd1(i, k) = 0.0
304 dnwd01(i, k) = 0.0
305 qcondc1(i, k) = 0.0
306 end do
307 end do
308
309 do i = 1, len
310 precip1(i) = 0.0
311 iflag1(i) = 0
312 wd1(i) = 0.0
313 cape1(i) = 0.0
314 VPrecip1(i, nd+1) = 0.0
315 end do
316
317 if (iflag_con.eq.3) then
318 do il = 1, len
319 sig1(il, nd) = sig1(il, nd) + 1.
320 sig1(il, nd) = min(sig1(il, nd), 12.1)
321 enddo
322 endif
323
324 !--------------------------------------------------------------------
325 ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
326 !--------------------------------------------------------------------
327
328 if (iflag_con.eq.3) then
329 CALL cv3_prelim(len, nd, nd + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, &
330 h1, hm1, th1)
331 endif
332
333 if (iflag_con.eq.4) then
334 CALL cv_prelim(len, nd, nd + 1, t1, q1, p1, ph1 &
335 , lv1, cpn1, tv1, gz1, h1, hm1)
336 endif
337
338 !--------------------------------------------------------------------
339 ! --- CONVECTIVE FEED
340 !--------------------------------------------------------------------
341
342 if (iflag_con.eq.3) then
343 CALL cv3_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1 &
344 , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! nd->na
345 endif
346
347 if (iflag_con.eq.4) then
348 CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1 &
349 , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)
350 endif
351
352 !--------------------------------------------------------------------
353 ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
354 ! (up through ICB for convect4, up through ICB+1 for convect3)
355 ! Calculates the lifted parcel virtual temperature at nk, the
356 ! actual temperature, and the adiabatic liquid water content.
357 !--------------------------------------------------------------------
358
359 if (iflag_con.eq.3) then
360 CALL cv3_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1 &
361 , tp1, tvp1, clw1, icbs1) ! nd->na
362 endif
363
364 if (iflag_con.eq.4) then
365 CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax &
366 , tp1, tvp1, clw1)
367 endif
368
369 !-------------------------------------------------------------------
370 ! --- TRIGGERING
371 !-------------------------------------------------------------------
372
373 if (iflag_con.eq.3) then
374 CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
375 buoybase1, iflag1, sig1, w01) ! nd->na
376 endif
377
378 if (iflag_con.eq.4) then
379 CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
380 endif
381
382 ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
383
384 ncum = 0
385 do i = 1, len
386 if(iflag1(i).eq.0)then
387 ncum = ncum+1
388 idcum(ncum) = i
389 endif
390 end do
391
392 IF (ncum.gt.0) THEN
393
394 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
395 ! --- COMPRESS THE FIELDS
396 ! (-> vectorization over convective gridpoints)
397 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
398
399 if (iflag_con.eq.3) then
400 CALL cv3_compress(len, nloc, ncum, nd, iflag1, nk1, icb1, icbs1, &
401 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
402 v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
403 sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
404 buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
405 tvp, clw, sig, w0)
406 endif
407
408 if (iflag_con.eq.4) then
409 CALL cv_compress( len, nloc, ncum, nd &
410 , iflag1, nk1, icb1 &
411 , cbmf1, plcl1, tnk1, qnk1, gznk1 &
412 , t1, q1, qs1, u1, v1, gz1 &
413 , h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1 &
414 , iflag, nk, icb &
415 , cbmf, plcl, tnk, qnk, gznk &
416 , t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw &
417 , dph )
418 endif
419
420 !-------------------------------------------------------------------
421 ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
422 ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
423 ! --- &
424 ! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
425 ! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
426 ! --- &
427 ! --- FIND THE LEVEL OF NEUTRAL BUOYANCY
428 !-------------------------------------------------------------------
429
430 if (iflag_con.eq.3) then
431 CALL cv3_undilute2(nloc, ncum, nd, icb, icbs, nk &
432 , tnk, qnk, gznk, t, q, qs, gz &
433 , p, h, tv, lv, pbase, buoybase, plcl &
434 , inb, tp, tvp, clw, hp, ep, sigp, buoy) !na->nd
435 endif
436
437 if (iflag_con.eq.4) then
438 CALL cv_undilute2(nloc, ncum, nd, icb, nk &
439 , tnk, qnk, gznk, t, q, qs, gz &
440 , p, dph, h, tv, lv &
441 , inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)
442 endif
443
444 !-------------------------------------------------------------------
445 ! --- CLOSURE
446 !-------------------------------------------------------------------
447
448 if (iflag_con.eq.3) then
449 CALL cv3_closure(nloc, ncum, nd, icb, inb &
450 , pbase, p, ph, tv, buoy &
451 , sig, w0, cape, m) ! na->nd
452 endif
453
454 if (iflag_con.eq.4) then
455 CALL cv_closure(nloc, ncum, nd, nk, icb &
456 , tv, tvp, p, ph, dph, plcl, cpn &
457 , iflag, cbmf)
458 endif
459
460 !-------------------------------------------------------------------
461 ! --- MIXING
462 !-------------------------------------------------------------------
463
464 if (iflag_con.eq.3) then
465 CALL cv3_mixing(nloc, ncum, nd, nd, icb, nk, inb, ph, t, q, &
466 qs, u, v, h, lv, qnk, hp, tv, tvp, ep, clw, m, sig, ment, &
467 qent, uent, vent, nent, sij, elij, ments, qents)
468 endif
469
470 if (iflag_con.eq.4) then
471 CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis &
472 , ph, t, q, qs, u, v, h, lv, qnk &
473 , hp, tv, tvp, ep, clw, cbmf &
474 , m, ment, qent, uent, vent, nent, sij, elij)
475 endif
476
477 !-------------------------------------------------------------------
478 ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
479 !-------------------------------------------------------------------
480
481 if (iflag_con.eq.3) then
482 CALL cv3_unsat(nloc, ncum, nd, nd, icb, inb &
483 , t, q, qs, gz, u, v, p, ph &
484 , th, tv, lv, cpn, ep, sigp, clw &
485 , m, ment, elij, delt, plcl &
486 , mp, qp, up, vp, wt, water, evap, b)! na->nd
487 endif
488
489 if (iflag_con.eq.4) then
490 CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph &
491 , h, lv, ep, sigp, clw, m, ment, elij &
492 , iflag, mp, qp, up, vp, wt, water, evap)
493 endif
494
495 !-------------------------------------------------------------------
496 ! --- YIELD
497 ! (tendencies, precipitation, variables of interface with other
498 ! processes, etc)
499 !-------------------------------------------------------------------
500
501 if (iflag_con.eq.3) then
502 CALL cv3_yield(nloc, ncum, nd, nd &
503 , icb, inb, delt &
504 , t, q, u, v, gz, p, ph, h, hp, lv, cpn, th &
505 , ep, clw, m, tp, mp, qp, up, vp &
506 , wt, water, evap, b &
507 , ment, qent, uent, vent, nent, elij, sig &
508 , tv, tvp &
509 , iflag, precip, VPrecip, ft, fq, fu, fv &
510 , upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc, wd)! na->nd
511 endif
512
513 if (iflag_con.eq.4) then
514 CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt &
515 , t, q, u, v, gz, p, ph, h, hp, lv, cpn &
516 , ep, clw, frac, m, mp, qp, up, vp &
517 , wt, water, evap &
518 , ment, qent, uent, vent, nent, elij &
519 , tv, tvp &
520 , iflag, wd, qprime, tprime &
521 , precip, cbmf, ft, fq, fu, fv, Ma, qcondc)
522 endif
523
524 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
525 ! --- passive tracers
526 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
527
528 if (iflag_con.eq.3) then
529 CALL cv3_tracer(nloc, len, ncum, nd, nd, &
530 ment, sij, da, phi)
531 endif
532
533 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
534 ! --- UNCOMPRESS THE FIELDS
535 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
536 ! set iflag1 = 42 for non convective points
537 do i = 1, len
538 iflag1(i) = 42
539 end do
540
541 if (iflag_con.eq.3) then
542 CALL cv3_uncompress(nloc, len, ncum, nd, idcum, iflag, precip, &
543 VPrecip, sig, w0, ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, &
544 qcondc, wd, cape, da, phi, mp, iflag1, precip1, VPrecip1, &
545 sig1, w01, ft1, fq1, fu1, fv1, inb1, Ma1, upwd1, dnwd1, &
546 dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
547 endif
548
549 if (iflag_con.eq.4) then
550 CALL cv_uncompress(nloc, len, ncum, nd, idcum &
551 , iflag &
552 , precip, cbmf &
553 , ft, fq, fu, fv &
554 , Ma, qcondc &
555 , iflag1 &
556 , precip1, cbmf1 &
557 , ft1, fq1, fu1, fv1 &
558 , Ma1, qcondc1 )
559 endif
560 ENDIF ! ncum>0
561
562 end SUBROUTINE cv_driver
563
564 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21