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

Last change on this file since 2225 was 2207, checked in by fhourdin, 10 years ago

Bi-gauussienne, suite.

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