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

Last change on this file since 2205 was 2205, checked in by fhourdin, 9 years ago

Introduction of a bi-gausian descrption of the subgrid scale
distribution of water associated with deep convection.

Introduction d'un schéma bi-gausien pour la distribution sous maille
de l'eau associée à la convection profonde.

Contrôlé par
iflag_cld_cv=2 (1 pour le schéma de Bony et Emanuel 2001)
autres paramètres de contrôle : tau_cld_cv et coefw_cld_cv

Arnaud Jam

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