/[lmdze]/trunk/libf/phylmd/cv_driver.f90
ViewVC logotype

Contents of /trunk/libf/phylmd/cv_driver.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (show annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 7 months ago) by guez
File size: 24012 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

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

  ViewVC Help
Powered by ViewVC 1.1.21