source: LMDZ5/trunk/libf/phylmd/cva_driver.F90 @ 2201

Last change on this file since 2201 was 2201, checked in by crio, 9 years ago

Expérience du tube de dentifrice: Ajout d'un terme de grande-échelle dans
la fermeture du schéma de convection. La fraction de la convergence
grande-échelle de masse au LFC ajoutée au flux de masse à la base calculé
par la fermeture en ALP est controlée par coef_clos_ls lu dans physiq.def
(valeur comprise entre 0 et 1, 0 par défaut).
Tube of toothpaste experiment: Additional large-scale term in the closure
formulation of the deep convection scheme. The fraction of the large-scale
convergence of mass at LFC additioned to the cloud-base mass-flux given by
the ALP closure is controlled by coef_clos_ls read in physiq.def (value
between 0 and 1, 0 by default).

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