source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phyiso/cv_driver.F90 @ 3488

Last change on this file since 3488 was 3331, checked in by acozic, 6 years ago

Add modification for isotopes

  • Property svn:executable set to *
File size: 44.8 KB
Line 
1
2! $Header$
3
4SUBROUTINE cv_driver(len, nd, ndp1, ntra, iflag_con, t1, q1, qs1, u1, v1, &
5    tra1, p1, ph1, iflag1, ft1, fq1, fu1, fv1, ftra1, precip1, vprecip1, &
6    cbmf1, sig1, w01, icb1, inb1, delt, ma1, upwd1, dnwd1, dnwd01, qcondc1, &
7    wd1, cape1, da1, phi1, mp1, phi21, d1a1, dam1, sij1, clw1, elij1, & !
8                                                                        ! RomP
9    evap1, ep1, epmlmmm1, eplamm1, & ! RomP
10    wdtraina1, wdtrainm1, & ! RomP
11    epmax_diag1  & ! epmax_cape
12#ifdef ISO
13     &                  ,xt1,fxt1,xtprecip1,xtvprecip1 &
14#ifdef DIAGISO
15     &          , water1,xtwater1,qp1,xtp1,xtevap1 &
16     &          , xtclw1 ! juste diagnostique
17     &          , wdtrain1,xtwdtrain1,tp1 &
18     &          , taux_cond_conv1 &
19     &          , fq_detrainement1,fq_ddft1,fq_fluxmasse1,fq_evapprecip1 &
20     &          , fxt_detrainement1,fxt_ddft1,fxt_fluxmasse1 &
21     &          , fxt_evapprecip1,Mi1,Amp_diag1,tg_save1 &
22     &          , f_detrainement1,q_detrainement1,xt_detrainement1 &
23#endif     
24#endif
25     &       )           
26
27  USE dimphy
28#ifdef ISO
29  USE infotrac_phy, ONLY: ntraciso,use_iso,niso,index_trac,ntraceurs_zone
30  USE isotopes_mod, ONLY: iso_eau,iso_HDO,ridicule,bidouille_anti_divergence
31#ifdef ISOVERIF
32    use isotopes_verif_mod, ONLY: errmax,errmaxrel,Tmin_verif,deltalim, &
33        iso_verif_egalite_choix,iso_verif_aberrant_choix, &
34        iso_verif_aberrant, iso_verif_positif,iso_verif_positif_nostop
35#endif
36#ifdef ISOTRAC
37    use isotrac_mod, only: option_traceurs,izone_ddft,izone_revap, &
38&        izone_poubelle, index_zone,option_tmin,izone_cond
39#ifdef ISOVERIF
40  USE isotopes_verif_mod, ONLY: iso_verif_traceur,iso_verif_traceur_justmass, &
41&       iso_verif_trac_masse_vect
42    use isotrac_routines_mod, ONLY: iso_verif_traceur_pbidouille
43#endif
44#endif
45#endif
46  IMPLICIT NONE
47
48  ! .............................START PROLOGUE............................
49
50
51  ! All argument names (except len,nd,ntra,nloc,delt and the flags) have a
52  ! "1" appended.
53  ! The "1" is removed for the corresponding compressed (local) variables.
54
55  ! PARAMETERS:
56  ! Name            Type         Usage            Description
57  ! ----------      ----------     -------  ----------------------------
58
59  ! len           Integer        Input        first (i) dimension
60  ! nd            Integer        Input        vertical (k) dimension
61  ! ndp1          Integer        Input        nd + 1
62  ! ntra          Integer        Input        number of tracors
63  ! iflag_con     Integer        Input        version of convect (3/4)
64  ! t1            Real           Input        temperature
65  ! q1            Real           Input        specific hum
66  ! qs1           Real           Input        sat specific hum
67  ! u1            Real           Input        u-wind
68  ! v1            Real           Input        v-wind
69  ! tra1          Real           Input        tracors
70  ! p1            Real           Input        full level pressure
71  ! ph1           Real           Input        half level pressure
72  ! iflag1        Integer        Output       flag for Emanuel conditions
73  ! ft1           Real           Output       temp tend
74  ! fq1           Real           Output       spec hum tend
75  ! fu1           Real           Output       u-wind tend
76  ! fv1           Real           Output       v-wind tend
77  ! ftra1         Real           Output       tracor tend
78  ! precip1       Real           Output       precipitation
79  ! VPrecip1      Real           Output       vertical profile of
80  ! precipitations
81  ! cbmf1         Real           Output       cloud base mass flux
82  ! sig1          Real           In/Out       section adiabatic updraft
83  ! w01           Real           In/Out       vertical velocity within adiab
84  ! updraft
85  ! delt          Real           Input        time step
86  ! Ma1           Real           Output       mass flux adiabatic updraft
87  ! upwd1         Real           Output       total upward mass flux
88  ! (adiab+mixed)
89  ! dnwd1         Real           Output       saturated downward mass flux
90  ! (mixed)
91  ! dnwd01        Real           Output       unsaturated downward mass flux
92  ! qcondc1       Real           Output       in-cld mixing ratio of
93  ! condensed water
94  ! wd1           Real           Output       downdraft velocity scale for
95  ! sfc fluxes
96  ! cape1         Real           Output       CAPE
97
98  ! wdtrainA1     Real           Output   precipitation detrained from
99  ! adiabatic draught;
100  ! used in tracer transport (cvltr)
101  ! wdtrainM1     Real           Output   precipitation detrained from mixed
102  ! draughts;
103  ! used in tracer transport (cvltr)
104  ! da1           Real           Output   used in tracer transport (cvltr)
105  ! phi1          Real           Output   used in tracer transport (cvltr)
106  ! mp1           Real           Output   used in tracer transport (cvltr)
107
108  ! phi21         Real           Output   used in tracer transport (cvltr)
109
110  ! d1a1          Real           Output   used in tracer transport (cvltr)
111  ! dam1          Real           Output   used in tracer transport (cvltr)
112
113  ! evap1         Real           Output
114  ! ep1           Real           Output
115  ! sij1        Real           Output
116  ! elij1         Real           Output
117
118  ! S. Bony, Mar 2002:
119  ! * Several modules corresponding to different physical processes
120  ! * Several versions of convect may be used:
121  ! - iflag_con=3: version lmd  (previously named convect3)
122  ! - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
123  ! + tard:     - iflag_con=5: version lmd with ice (previously named convectg)
124  ! S. Bony, Oct 2002:
125  ! * Vectorization of convect3 (ie version lmd)
126
127  ! ..............................END PROLOGUE.............................
128
129
130  ! Input
131  INTEGER len
132  INTEGER nd
133  INTEGER ndp1
134  INTEGER noff
135  INTEGER iflag_con
136  INTEGER ntra
137  REAL delt
138  REAL t1(len, nd)
139  REAL q1(len, nd)
140  REAL qs1(len, nd)
141  REAL u1(len, nd)
142  REAL v1(len, nd)
143  REAL tra1(len, nd, ntra)
144  REAL p1(len, nd)
145  REAL ph1(len, ndp1)
146
147  ! Output
148  INTEGER iflag1(len)
149  REAL ft1(len, nd)
150  REAL fq1(len, nd)
151  REAL fu1(len, nd)
152  REAL fv1(len, nd)
153  REAL ftra1(len, nd, ntra)
154  REAL precip1(len)
155  REAL cbmf1(len)
156  REAL sig1(klon, klev)
157  REAL w01(klon, klev)
158  REAL vprecip1(len, nd+1)
159  REAL evap1(len, nd) !RomP
160  REAL ep1(len, nd) !RomP
161  REAL ma1(len, nd)
162  REAL upwd1(len, nd)
163  REAL dnwd1(len, nd)
164  REAL dnwd01(len, nd)
165
166  REAL qcondc1(len, nd) ! cld
167  REAL wd1(len) ! gust
168  REAL cape1(len)
169
170  ! RomP >>>
171  REAL wdtraina1(len, nd), wdtrainm1(len, nd)
172  REAL sij1(len, nd, nd), elij1(len, nd, nd)
173  REAL da1(len, nd), phi1(len, nd, nd), mp1(len, nd)
174
175  REAL phi21(len, nd, nd)
176  REAL d1a1(len, nd), dam1(len, nd)
177  REAL epmlmmm1(len, nd, nd), eplamm1(len, nd)
178  ! RomP <<<
179  REAL epmax_diag1 (len) ! epmax_cape     
180
181
182#ifdef ISO
183      !integer niso
184      real xt1(ntraciso,len,nd)
185      real fxt1(ntraciso,len,nd)
186      real xtprecip1(ntraciso,len)
187      real xtvprecip1(ntraciso,len,nd+1)
188#ifdef DIAGISO
189      real water1(len,nd)
190      real xtwater1(niso,len,nd)
191      real qp1(len,nd),xtp1(ntraciso,len,nd)
192      real xtevap1(niso,len,nd)
193      real wdtrain1(len,nd)
194      real xtwdtrain1(niso,len,nd)
195      real taux_cond_conv1(len,nd)
196      real epmax_diag1(len)
197      real fq_detrainement1(len,nd)
198      real f_detrainement1(len,nd)
199      real q_detrainement1(len,nd)
200      real fq_ddft1(len,nd)
201      real fq_fluxmasse1(len,nd)
202      real Amp_diag1(len,nd)
203!      real tg_save1(len,nd)
204      real fq_evapprecip1(len,nd)
205      real fxt_detrainement1(niso,len,nd)
206      real xt_detrainement1(niso,len,nd)
207      real fxt_ddft1(niso,len,nd)
208      real fxt_fluxmasse1(niso,len,nd)
209      real fxt_evapprecip1(niso,len,nd)
210      real Mi1(len,nd)
211!      real mentbas1(len,nd)
212!      real qentbas1(len,nd), xtentbas1(niso,len,nd)
213#endif
214#ifdef ISOVERIF     
215      real deltaD
216#endif     
217#endif
218
219  ! -------------------------------------------------------------------
220  ! Original Prologue by Kerry Emanuel.
221  ! -------------------------------------------------------------------
222  ! --- ARGUMENTS
223  ! -------------------------------------------------------------------
224  ! --- On input:
225
226  ! t:   Array of absolute temperature (K) of dimension ND, with first
227  ! index corresponding to lowest model level. Note that this array
228  ! will be altered by the subroutine if dry convective adjustment
229  ! occurs and if IPBL is not equal to 0.
230
231  ! q:   Array of specific humidity (gm/gm) of dimension ND, with first
232  ! index corresponding to lowest model level. Must be defined
233  ! at same grid levels as T. Note that this array will be altered
234  ! if dry convective adjustment occurs and if IPBL is not equal to 0.
235
236  ! qs:  Array of saturation specific humidity of dimension ND, with first
237  ! index corresponding to lowest model level. Must be defined
238  ! at same grid levels as T. Note that this array will be altered
239  ! if dry convective adjustment occurs and if IPBL is not equal to 0.
240
241  ! u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
242  ! index corresponding with the lowest model level. Defined at
243  ! same levels as T. Note that this array will be altered if
244  ! dry convective adjustment occurs and if IPBL is not equal to 0.
245
246  ! v:   Same as u but for meridional velocity.
247
248  ! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
249  ! where NTRA is the number of different tracers. If no
250  ! convective tracer transport is needed, define a dummy
251  ! input array of dimension (ND,1). Tracers are defined at
252  ! same vertical levels as T. Note that this array will be altered
253  ! if dry convective adjustment occurs and if IPBL is not equal to 0.
254
255  ! p:   Array of pressure (mb) of dimension ND, with first
256  ! index corresponding to lowest model level. Must be defined
257  ! at same grid levels as T.
258
259  ! ph:  Array of pressure (mb) of dimension ND+1, with first index
260  ! corresponding to lowest level. These pressures are defined at
261  ! levels intermediate between those of P, T, Q and QS. The first
262  ! value of PH should be greater than (i.e. at a lower level than)
263  ! the first value of the array P.
264
265  ! nl:  The maximum number of levels to which convection can penetrate, plus
266  ! 1.
267  ! NL MUST be less than or equal to ND-1.
268
269  ! delt: The model time step (sec) between calls to CONVECT
270
271  ! ----------------------------------------------------------------------------
272  ! ---   On Output:
273
274  ! iflag: An output integer whose value denotes the following:
275  ! VALUE   INTERPRETATION
276  ! -----   --------------
277  ! 0     Moist convection occurs.
278  ! 1     Moist convection occurs, but a CFL condition
279  ! on the subsidence warming is violated. This
280  ! does not cause the scheme to terminate.
281  ! 2     Moist convection, but no precip because ep(inb) lt 0.0001
282  ! 3     No moist convection because new cbmf is 0 and old cbmf is 0.
283  ! 4     No moist convection; atmosphere is not
284  ! unstable
285  ! 6     No moist convection because ihmin le minorig.
286  ! 7     No moist convection because unreasonable
287  ! parcel level temperature or specific humidity.
288  ! 8     No moist convection: lifted condensation
289  ! level is above the 200 mb level.
290  ! 9     No moist convection: cloud base is higher
291  ! then the level NL-1.
292
293  ! ft:   Array of temperature tendency (K/s) of dimension ND, defined at
294  ! same
295  ! grid levels as T, Q, QS and P.
296
297  ! fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
298  ! defined at same grid levels as T, Q, QS and P.
299
300  ! fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
301  ! defined at same grid levels as T.
302
303  ! fv:   Same as FU, but for forcing of meridional velocity.
304
305  ! ftra: Array of forcing of tracer content, in tracer mixing ratio per
306  ! second, defined at same levels as T. Dimensioned (ND,NTRA).
307
308  ! precip: Scalar convective precipitation rate (mm/day).
309
310  ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).
311
312  ! wd:   A convective downdraft velocity scale. For use in surface
313  ! flux parameterizations. See convect.ps file for details.
314
315  ! tprime: A convective downdraft temperature perturbation scale (K).
316  ! For use in surface flux parameterizations. See convect.ps
317  ! file for details.
318
319  ! qprime: A convective downdraft specific humidity
320  ! perturbation scale (gm/gm).
321  ! For use in surface flux parameterizations. See convect.ps
322  ! file for details.
323
324  ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
325  ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
326  ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"
327  ! by the calling program between calls to CONVECT.
328
329  ! det:   Array of detrainment mass flux of dimension ND.
330
331  ! -------------------------------------------------------------------
332
333  ! Local arrays
334
335
336  INTEGER i, k, n, il, j
337  INTEGER icbmax
338  INTEGER nk1(klon)
339  INTEGER icb1(klon)
340  INTEGER inb1(klon)
341  INTEGER icbs1(klon)
342
343  REAL plcl1(klon)
344  REAL tnk1(klon)
345  REAL qnk1(klon)
346  REAL gznk1(klon)
347  REAL pnk1(klon)
348  REAL qsnk1(klon)
349  REAL pbase1(klon)
350  REAL buoybase1(klon)
351
352  REAL lv1(klon, klev)
353  REAL cpn1(klon, klev)
354  REAL tv1(klon, klev)
355  REAL gz1(klon, klev)
356  REAL hm1(klon, klev)
357  REAL h1(klon, klev)
358  REAL tp1(klon, klev)
359  REAL tvp1(klon, klev)
360  REAL clw1(klon, klev)
361  REAL th1(klon, klev)
362
363#ifdef ISO
364      integer ixt
365      real xtclw1(ntraciso,klon,klev)
366      real tg_save1(klon,klev)
367      REAL xtnk1(ntraciso,klon)
368#endif
369
370  INTEGER ncum
371
372  ! (local) compressed fields:
373
374  ! ym      integer nloc
375  ! ym      parameter (nloc=klon) ! pour l'instant
376#define nloc klon
377  INTEGER idcum(nloc)
378  INTEGER iflag(nloc), nk(nloc), icb(nloc)
379  INTEGER nent(nloc, klev)
380  INTEGER icbs(nloc)
381  INTEGER inb(nloc), inbis(nloc)
382
383  REAL cbmf(nloc), plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)
384  REAL t(nloc, klev), q(nloc, klev), qs(nloc, klev)
385  REAL u(nloc, klev), v(nloc, klev)
386  REAL gz(nloc, klev), h(nloc, klev), lv(nloc, klev), cpn(nloc, klev)
387  REAL p(nloc, klev), ph(nloc, klev+1), tv(nloc, klev), tp(nloc, klev)
388  REAL clw(nloc, klev)
389  REAL dph(nloc, klev)
390  REAL pbase(nloc), buoybase(nloc), th(nloc, klev)
391  REAL tvp(nloc, klev)
392  REAL sig(nloc, klev), w0(nloc, klev)
393  REAL hp(nloc, klev), ep(nloc, klev), sigp(nloc, klev)
394  REAL frac(nloc), buoy(nloc, klev)
395  REAL cape(nloc)
396  REAL m(nloc, klev), ment(nloc, klev, klev), qent(nloc, klev, klev)
397  REAL uent(nloc, klev, klev), vent(nloc, klev, klev)
398  REAL ments(nloc, klev, klev), qents(nloc, klev, klev)
399  REAL sij(nloc, klev, klev), elij(nloc, klev, klev)
400  REAL qp(nloc, klev), up(nloc, klev), vp(nloc, klev)
401  REAL wt(nloc, klev), water(nloc, klev), evap(nloc, klev)
402  REAL b(nloc, klev), ft(nloc, klev), fq(nloc, klev)
403  REAL fu(nloc, klev), fv(nloc, klev)
404  REAL upwd(nloc, klev), dnwd(nloc, klev), dnwd0(nloc, klev)
405  REAL ma(nloc, klev), mike(nloc, klev), tls(nloc, klev)
406  REAL tps(nloc, klev), qprime(nloc), tprime(nloc)
407  REAL precip(nloc)
408  REAL vprecip(nloc, klev+1)
409  REAL tra(nloc, klev, ntra), trap(nloc, klev, ntra)
410  REAL ftra(nloc, klev, ntra), traent(nloc, klev, klev, ntra)
411  REAL qcondc(nloc, klev) ! cld
412  REAL wd(nloc) ! gust
413
414  ! RomP >>>
415  REAL da(nloc, klev), phi(nloc, klev, klev), mp(nloc, klev)
416  REAL epmlmmm(nloc, klev, klev), eplamm(nloc, klev)
417  REAL phi2(nloc, klev, klev)
418  REAL d1a(nloc, klev), dam(nloc, klev)
419  REAL wdtraina(nloc, klev), wdtrainm(nloc, klev)
420  REAL sigd(nloc)
421  ! RomP <<<
422  REAL epmax_diag(nloc) ! epmax_cape
423#ifdef ISO
424      real xt(ntraciso,nloc,klev)
425      real xtnk(ntraciso,nloc)
426      real xtclw(ntraciso,nloc,klev)
427      real xtp(ntraciso,nloc,klev)
428      real xtent(ntraciso,nloc,klev,klev)
429      real xtelij(ntraciso,nloc,klev,klev)
430!      real xtep(ntraciso,nloc,klev) ! le 7 mai: on supprime xtep
431      real xtwater(ntraciso,nloc,klev)
432      real xtevap(ntraciso,nloc,klev)
433      real fxt(ntraciso,nloc,klev)
434      real xtprecip(ntraciso,nloc)
435      real tg_save(nloc,klev) ! températures de condensation
436      real xtvprecip(ntraciso,nloc,klev+1)
437#ifdef DIAGISO
438      real wdtrain(nloc,klev)
439      real xtwdtrain(niso,nloc,klev)
440      real taux_cond_conv(nloc,klev)
441      real epmax_diag(nloc)
442      real fxt_detrainement(niso,nloc,klev)
443      real fxt_fluxmasse(niso,nloc,klev)
444      real fxt_evapprecip(niso,nloc,klev)
445      real fxt_ddft(niso,nloc,klev)
446      real fq_detrainement(nloc,klev)
447      real f_detrainement(nloc,klev)
448      real q_detrainement(nloc,klev)
449      real xt_detrainement(niso,nloc,klev)
450      real fq_fluxmasse(nloc,klev)
451      real Amp_diag(nloc,klev)
452      real fq_evapprecip(nloc,klev)
453      real fq_ddft(nloc,klev)
454!      real Mi(nloc,klev) déjà sorti dans cv3_yield: m
455!      real mentbas(nloc,klev): déjà sortie: ment
456!      real qentbas(nloc,klev): déjà sortie: qent
457#endif
458#ifdef ISOTRAC
459      integer iiso,ixt_ddft,ixt_poubelle,ixt_revap
460      integer izone
461!#ifdef ISOVERIF
462!      integer iso_verif_positif_nostop
463!#endif     
464#endif
465#endif
466
467
468#ifdef ISO
469        ! verif     
470#ifdef ISOVERIF
471       write(*,*) 'cv_driver 391: ENTREE dans cv_driver'
472!        write(*,*) 'nloc=',nloc
473       if (use_iso(iso_eau)) then
474        do k = 1, klev
475         do i = 1, nloc         
476            call iso_verif_egalite_choix(xt1(iso_eau,i,k),q1(i,k), &
477     &           'cv_driver 343',errmax,errmaxrel)           
478         enddo
479        enddo
480       endif !if (use_iso(iso_eau)) then
481       if (use_iso(iso_HDO)) then
482        do k = 1, klev
483         do i = 1, nloc 
484           if (q1(i,k).gt.ridicule) then       
485            call iso_verif_aberrant(xt1(iso_hdo,i,k)/q1(i,k), &
486     &           'cv_driver 362')     
487           endif !if (q1(i,k).gt.ridicule) then       
488         enddo
489        enddo
490       endif !if (use_iso(iso_eau)) then
491#ifdef ISOTRAC
492        do k = 1, klev
493        do i = 1, klon 
494           call iso_verif_traceur(xt1(1,i,k),'cv_driver 413')
495        enddo
496        enddo
497#endif 
498#endif
499       ! end verif       
500#endif
501
502
503  nent(:, :) = 0
504  ! -------------------------------------------------------------------
505  ! --- SET CONSTANTS AND PARAMETERS
506  ! -------------------------------------------------------------------
507  ! print *, '-> cv_driver'      !jyg
508  ! -- set simulation flags:
509  ! (common cvflag)
510
511  CALL cv_flag(0)
512
513  ! -- set thermodynamical constants:
514  ! (common cvthermo)
515
516  CALL cv_thermo(iflag_con)
517
518  ! -- set convect parameters
519
520  ! includes microphysical parameters and parameters that
521  ! control the rate of approach to quasi-equilibrium)
522  ! (common cvparam)
523
524
525  IF (iflag_con==30) THEN
526    CALL cv30_param(nd, delt)
527  END IF
528
529  IF (iflag_con==4) THEN
530    CALL cv_param(nd)
531#ifdef ISO
532       write(*,*) 'cv_driver 454: isos pas prévus ici'
533       stop
534#endif
535  END IF
536
537  ! ---------------------------------------------------------------------
538  ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
539  ! ---------------------------------------------------------------------
540
541  ft1(:, :) = 0.0
542  fq1(:, :) = 0.0
543  fu1(:, :) = 0.0
544  fv1(:, :) = 0.0
545  tvp1(:, :) = 0.0
546  tp1(:, :) = 0.0
547  clw1(:, :) = 0.0
548  ! ym
549  clw(:, :) = 0.0
550  gz1(:, :) = 0.
551  vprecip1(:, :) = 0.
552  ma1(:, :) = 0.0
553  upwd1(:, :) = 0.0
554  dnwd1(:, :) = 0.0
555  dnwd01(:, :) = 0.0
556  qcondc1(:, :) = 0.0
557
558  ftra1(:, :, :) = 0.0
559
560  elij1(:, :, :) = 0.0
561  sij1(:, :, :) = 0.0
562
563  precip1(:) = 0.0
564  iflag1(:) = 0
565  wd1(:) = 0.0
566  cape1(:) = 0.0
567  epmax_diag1(:) = 0.0 ! epmax_cape
568
569
570#ifdef ISOVERIF
571#ifdef ISOTRAC     
572        do k=1,klev       
573         do i=1,klon       
574           call iso_verif_traceur(xt1(1,i,k),'cv_driver 510')
575         enddo
576       enddo
577#endif
578#endif         
579
580#ifdef ISO
581      fxt1(:,:,:)=0.0
582      xtclw1(:,:,:)=0.0
583      xtclw(:,:,:)=0.0
584      xtvprecip1(:,:,:) = 0.
585      tg_save1(:,:) = 0.
586#ifdef DIAGISO 
587      do k=1,nd
588       do i=1,len
589        do ixt=1,niso       
590         xtwater1(ixt,i,k)=0.0
591         xtwdtrain1(ixt,i,k)=0.0
592         fxt_detrainement1(ixt,i,k)=0.0
593         xt_detrainement1(ixt,i,k)=0.0
594         fxt_ddft1(ixt,i,k)=0.0
595         fxt_fluxmasse1(ixt,i,k)=0.0
596         fxt_evapprecip1(ixt,i,k)=0.0
597         xtp1(ixt,i,k)=0.0
598         xtevap1(ixt,i,k)=0.0
599!         xtentbas(ixt,i,k)=0.0
600        enddo
601#ifdef ISOTRAC
602        do ixt=1+niso,ntraciso
603           xtp1(ixt,i,k)=0.0
604        enddo !do ixt=1+niso,ntraciso
605#endif       
606        evap1(i,k)=0.0
607        wdtrain1(i,k)=0.0
608        taux_cond_conv1(i,k)=0.0
609        water1(i,k)=0.0
610        Mp1(i,k)=0.0
611        qp1(i,k)=0.0
612        fq_detrainement1(i,k)=0.0
613        fq_ddft1(i,k)=0.0
614        fq_fluxmasse1(i,k)=0.0
615        Amp_diag1(i,k)=0.0
616        tg_save1(i,k)=t1(i,k)
617        fq_evapprecip1(i,k)=0.0
618        Mi1(i,k)=0.0
619!        qentbas1(i,k)=0.0
620!        mentbas1(i,k)=0.0
621        enddo !do i=1,len
622       enddo !do k=1,nd
623#endif         
624#endif
625
626  IF (iflag_con==30) THEN
627    DO il = 1, len
628      sig1(il, nd) = sig1(il, nd) + 1.
629      sig1(il, nd) = amin1(sig1(il,nd), 12.1)
630    END DO
631  END IF
632
633  ! RomP >>>
634  wdtraina1(:, :) = 0.
635  wdtrainm1(:, :) = 0.
636  da1(:, :) = 0.
637  phi1(:, :, :) = 0.
638  epmlmmm1(:, :, :) = 0.
639  eplamm1(:, :) = 0.
640  mp1(:, :) = 0.
641  evap1(:, :) = 0.
642  ep1(:, :) = 0.
643  sij1(:, :, :) = 0.
644  elij1(:, :, :) = 0.
645  phi21(:, :, :) = 0.
646  d1a1(:, :) = 0.
647  dam1(:, :) = 0.
648  ! RomP <<<
649
650  ! --------------------------------------------------------------------
651  ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
652  ! --------------------------------------------------------------------
653
654  IF (iflag_con==30) THEN
655
656    ! print*,'Emanuel version 30 '
657    CALL cv30_prelim(len, nd, ndp1, t1, q1, p1, ph1 & ! nd->na
658      , lv1, cpn1, tv1, gz1, h1, hm1, th1)
659  END IF
660
661  IF (iflag_con==4) THEN
662    CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, &
663      hm1)
664  END IF
665
666  ! --------------------------------------------------------------------
667  ! --- CONVECTIVE FEED
668  ! --------------------------------------------------------------------
669
670  IF (iflag_con==30) THEN
671    CALL cv30_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1 & !
672                                                             ! nd->na
673      , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1 &
674#ifdef ISO
675      ,xt1,xtnk1 &   
676#endif
677    )
678  END IF
679
680  IF (iflag_con==4) THEN
681    CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1, nk1, icb1, icbmax, &
682      iflag1, tnk1, qnk1, gznk1, plcl1)
683  END IF
684
685  ! --------------------------------------------------------------------
686  ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
687  ! (up through ICB for convect4, up through ICB+1 for convect3)
688  ! Calculates the lifted parcel virtual temperature at nk, the
689  ! actual temperature, and the adiabatic liquid water content.
690  ! --------------------------------------------------------------------
691
692  IF (iflag_con==30) THEN
693
694#ifdef ISOVERIF
695       write(*,*) 'cv_driver 593: avant cv3_undilute1'
696       do k=1,klev       
697         do i=1,klon
698          if (use_iso(iso_HDO)) then
699           if (q1(i,k).gt.ridicule) then
700            call iso_verif_aberrant(xt1(iso_hdo,i,k)/q1(i,k), &
701     &            'cv_driver 502')
702            endif ! if (q1(i,k).gt.ridicule) then
703            endif !if (use_iso(iso_HDO)) then
704#ifdef ISOTRAC
705           call iso_verif_traceur(xt1(1,i,k),'cv_driver 601')
706#endif     
707         enddo
708       enddo       
709#endif
710
711
712    CALL cv30_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1 & ! nd->na
713      , tp1, tvp1, clw1, icbs1 &
714#ifdef ISO
715     &        ,xt1,xtclw1,tg_save1 &
716#endif
717     &        )
718
719
720#ifdef ISO
721!c--debug
722#ifdef ISOVERIF
723       write(*,*) 'cv_driver 621: après cv3_undilute1'
724       do k = 1, klev
725        do i = 1, klon
726         if (use_iso(iso_eau)) then
727         call iso_verif_egalite_choix(xtclw1(iso_eau,i,k),clw1(i,k), &
728     &           'undil1',errmax,errmaxrel)
729         call iso_verif_egalite_choix(xt1(iso_eau,i,k),q1(i,k), &
730     &           'undil1',errmax,errmaxrel)
731         endif !if (use_iso(iso_eau)) then
732#ifdef ISOTRAC
733           call iso_verif_traceur(xt1(1,i,k),'cv_driver 623')
734#endif           
735        enddo
736       enddo
737#endif
738       !write(*,*) 'SORTIE DE CV3_UNDILUTE1'
739#endif
740
741  END IF
742
743  IF (iflag_con==4) THEN
744    CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, tp1, &
745      tvp1, clw1)
746  END IF
747
748  ! -------------------------------------------------------------------
749  ! --- TRIGGERING
750  ! -------------------------------------------------------------------
751
752  IF (iflag_con==30) THEN
753    CALL cv30_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1 & !
754                                                                 ! nd->na
755      , pbase1, buoybase1, iflag1, sig1, w01)
756  END IF
757
758  IF (iflag_con==4) THEN
759    CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
760  END IF
761
762  ! =====================================================================
763  ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
764  ! =====================================================================
765
766  ncum = 0
767  DO i = 1, len
768    IF (iflag1(i)==0) THEN
769      ncum = ncum + 1
770      idcum(ncum) = i
771    END IF
772  END DO
773
774  ! print*,'cv_driver : klon, ncum = ',len,ncum
775
776  IF (ncum>0) THEN
777
778    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
779    ! --- COMPRESS THE FIELDS
780    ! (-> vectorization over convective gridpoints)
781    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
782
783    IF (iflag_con==30) THEN
784#ifdef ISO
785#ifdef ISOVERIF
786       write(*,*) 'cv_driver 771: tg_save1(1,1)=',tg_save1(1,1)       
787        write(*,*) 'xt1(iso_eau,1,1),q1(1,1)=',xt1(iso_eau,1,1),q1(1,1)
788        !write(*,*) 'xt1(iso_eau,14,1),q1(14,1)=',xt1(iso_eau,14,1),q1(14,1)
789        write(*,*) 'iso_eau,use_iso=',iso_eau,use_iso
790       do k = 1, klev
791        do i = 1, nloc
792        if (use_iso(iso_eau)) then
793            call iso_verif_egalite_choix(xtclw1(iso_eau,i,k),clw1(i,k), &
794     &            'cv_driver 541a',errmax,errmaxrel)
795            call iso_verif_egalite_choix(xt1(iso_eau,i,k),q1(i,k), &
796     &            'cv_driver 541b',errmax,errmaxrel)
797        endif !  if (use_iso(iso_eau)) then
798#ifdef ISOTRAC
799           call iso_verif_traceur(xt1(1,i,k),'cv_driver 689')
800#endif             
801        enddo
802       enddo   
803#endif
804#endif
805
806      CALL cv30_compress(len, nloc, ncum, nd, ntra, iflag1, nk1, icb1, icbs1, &
807        plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, &
808        gz1, th1, tra1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, &
809        w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, &
810        q, qs, u, v, gz, th, tra, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, &
811        w0 &
812#ifdef ISO
813     &    ,xtnk1,xt1,xtclw1,tg_save1 &
814     &    ,xtnk,xt,xtclw,tg_save &
815#endif
816     &    )
817
818#ifdef ISO
819#ifdef ISOVERIF
820       write(*,*) 'cv_driver 720: après cv3_compress'           
821       write(*,*) 'tg_save(1,1)=',tg_save(1,1) 
822       do k = 1, klev
823        do i = 1, ncum
824         if (use_iso(iso_eau)) then
825            call iso_verif_egalite_choix(xtclw(iso_eau,i,k),clw(i,k), &
826     &            'cv_driver 598',errmax,errmaxrel)
827            call iso_verif_egalite_choix(xt(iso_eau,i,k),q(i,k), &
828     &            'cv_driver 600',errmax,errmaxrel)
829         endif !  if (use_iso(iso_eau)) then
830         if (use_iso(iso_HDO)) then   
831              call iso_verif_aberrant_choix( &
832     &            xt(iso_HDO,i,k),q(i,k), &
833     &            ridicule,deltalim,'cv_driver 735, apres compress')
834         endif !if (use_iso(iso_HDO)) then
835#ifdef ISOTRAC
836           call iso_verif_traceur(xt(1,i,k),'cv_driver 726')
837#endif                   
838         call iso_verif_positif(tg_save(i,k)-Tmin_verif, &
839     &        'cv_driver 806')   
840        enddo
841       enddo
842       if (use_iso(iso_eau)) then
843       do i = 1, ncum
844            call iso_verif_egalite_choix(xtnk(iso_eau,i),qnk(i), &
845     &            'cv_driver 834',errmax,errmaxrel)
846       enddo
847       endif !if (use_iso(iso_eau)) then
848#endif
849#endif
850
851    END IF
852
853    IF (iflag_con==4) THEN
854      CALL cv_compress(len, nloc, ncum, nd, iflag1, nk1, icb1, cbmf1, plcl1, &
855        tnk1, qnk1, gznk1, t1, q1, qs1, u1, v1, gz1, h1, lv1, cpn1, p1, ph1, &
856        tv1, tp1, tvp1, clw1, iflag, nk, icb, cbmf, plcl, tnk, qnk, gznk, t, &
857        q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw, dph)
858    END IF
859
860    ! -------------------------------------------------------------------
861    ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
862    ! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
863    ! ---   &
864    ! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
865    ! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
866    ! ---   &
867    ! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
868    ! -------------------------------------------------------------------
869
870    IF (iflag_con==30) THEN
871      CALL cv30_undilute2(nloc, ncum, nd, icb, icbs, nk & !na->nd
872        , tnk, qnk, gznk, t, q, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, &
873        inb, tp, tvp, clw, hp, ep, sigp, buoy &
874#ifdef ISO
875     &  ,xtnk,xt,xtclw,tg_save &
876#endif
877     &   )
878
879#ifdef ISO
880#ifdef ISOVERIF
881        write(*,*) 'cv_driver 863: apres call cv30_undilute2'
882        write(*,*) 'tg_save(1,1)=',tg_save(1,1)
883       do k = 1, klev
884        do i = 1, ncum
885         if (use_iso(iso_eau)) then
886            call iso_verif_egalite_choix(xtclw(iso_eau,i,k),clw(i,k), &
887     &            'cv_driver 650',errmax,errmaxrel)
888            call iso_verif_egalite_choix(xt(iso_eau,i,k),q(i,k), &
889     &            'cv_driver 652',errmax,errmaxrel)
890         endif !  if (use_iso(iso_eau)) then
891         if (use_iso(iso_HDO)) then   
892              call iso_verif_aberrant_choix( &
893     &            xt(iso_HDO,i,k),q(i,k), &
894     &            ridicule,deltalim,'cv_driver 794, apres undilute2')
895         endif !if (use_iso(iso_HDO)) then 
896         call iso_verif_positif(tg_save(i,k)-Tmin_verif, &
897     &        'cv_driver 859')
898#ifdef ISOTRAC
899           call iso_verif_traceur(xt(1,i,k),'cv_driver 780')
900           call iso_verif_traceur(xtclw(1,i,k),'cv_driver 781')
901#endif               
902        enddo
903       enddo
904#endif
905       !write(*,*) 'SORTIE CV3_UNDILUTE2'
906#endif
907    END IF
908
909    IF (iflag_con==4) THEN
910      CALL cv_undilute2(nloc, ncum, nd, icb, nk, tnk, qnk, gznk, t, q, qs, &
911        gz, p, dph, h, tv, lv, inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)
912    END IF
913
914    ! -------------------------------------------------------------------
915    ! --- CLOSURE
916    ! -------------------------------------------------------------------
917
918    IF (iflag_con==30) THEN
919      CALL cv30_closure(nloc, ncum, nd, icb, inb & ! na->nd
920        , pbase, p, ph, tv, buoy, sig, w0, cape, m)
921
922      ! epmax_cape
923      call cv30_epmax_fn_cape(nloc,ncum,nd &
924                ,cape,ep,hp,icb,inb,clw,nk,t,h,lv &
925                ,epmax_diag)
926        ! on écrase ep et recalcule hp
927    END IF
928
929    IF (iflag_con==4) THEN
930      CALL cv_closure(nloc, ncum, nd, nk, icb, tv, tvp, p, ph, dph, plcl, &
931        cpn, iflag, cbmf)
932    END IF
933   
934
935    ! -------------------------------------------------------------------
936    ! --- MIXING
937    ! -------------------------------------------------------------------
938
939    IF (iflag_con==30) THEN
940      CALL cv30_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb & !
941                                                                ! na->nd
942        , ph, t, q, qs, u, v, tra, h, lv, qnk, hp, tv, tvp, ep, clw, m, sig, &
943        ment, qent, uent, vent, sij, elij, ments, qents, traent &
944#ifdef ISO
945     &                     ,xt,xtnk,xtclw &
946     &                     ,xtent,xtelij,tg_save &
947#endif
948     &    )
949
950
951#ifdef ISO
952#ifdef ISOVERIF
953       write(*,*) 'cv_driver 837: après cv3_mixing'
954       do k = 1, klev
955       do j = 1, klev
956        do i = 1, ncum
957         if (use_iso(iso_eau)) then
958            call iso_verif_egalite_choix(xtelij(iso_eau,i,j,k), &
959     &            elij(i,j,k),'cv_driver 881',errmax,errmaxrel)
960            call iso_verif_egalite_choix(xtent(iso_eau,i,j,k), &
961     &            qent(i,j,k),'cv_driver 882',errmax,errmaxrel)
962         endif !  if (use_iso(iso_eau)) then
963#ifdef ISOTRAC
964           call iso_verif_traceur_justmass(xtent(1,i,j,k), &
965     &           'cv_driver 846')
966           call iso_verif_traceur_justmass(xtelij(1,i,j,k), &
967     &           'cv_driver 847')
968           ! on ne vérfier pas le deltaD ici car peut dépasser le seuil
969           ! raisonable pour températures très froides.
970#endif               
971        enddo
972       enddo
973       enddo !do k = 1, klev
974       do k = 1, klev
975        do i = 1, ncum
976         if (use_iso(iso_eau)) then
977            call iso_verif_egalite_choix(xt(iso_eau,i,k),q(i,k), &
978     &            'cv_driver 709',errmax,errmaxrel)
979         endif !  if (use_iso(iso_eau)) then
980#ifdef ISOTRAC
981           call iso_verif_traceur(xt(1,i,k),'cv_driver 856')
982           if (option_tmin.eq.1) then
983             if (iso_verif_positif_nostop(xtclw(index_trac( &
984     &           izone_cond,iso_eau),i,k)-xtclw(iso_eau,i,k) &
985     &           ,'cv_driver 909').eq.1) then
986               write(*,*) 'i,k=',i,k
987               write(*,*) 'xtclw=',xtclw(:,i,k)
988               stop
989             endif !if (iso_verif_positif_nostop(xtclw(index_trac(
990           endif !if ((option_traceurs.eq.17).or.
991#endif 
992#ifdef DIAGISO
993        call iso_verif_positif(tg_save(i,k)-50.0,'cv_driver 926')
994#endif
995        enddo
996       enddo !do k = 1, klev 
997   
998#endif
999#endif
1000
1001    END IF
1002
1003    IF (iflag_con==4) THEN
1004      CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis, ph, t, q, qs, u, v, &
1005        h, lv, qnk, hp, tv, tvp, ep, clw, cbmf, m, ment, qent, uent, vent, &
1006        nent, sij, elij)
1007    END IF
1008
1009    ! -------------------------------------------------------------------
1010    ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
1011    ! -------------------------------------------------------------------
1012
1013    IF (iflag_con==30) THEN
1014
1015#ifdef ISO
1016#ifdef ISOVERIF
1017       do k = 1, klev
1018        do i = 1, ncum
1019         if (use_iso(iso_eau)) then
1020            call iso_verif_egalite_choix(xt(iso_eau,i,k),q(i,k), &
1021     &            'cv_driver 753',errmax,errmaxrel)
1022         endif !  if (use_iso(iso_eau)) then
1023#ifdef ISOTRAC
1024           call iso_verif_traceur(xt(1,i,k),'cv_driver 885')           
1025#endif               
1026        enddo
1027       enddo !do k = 1, klev     
1028#endif
1029#endif
1030
1031      ! RomP >>>
1032      CALL cv30_unsat(nloc, ncum, nd, nd, ntra, icb, inb & ! na->nd
1033        , t, q, qs, gz, u, v, tra, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, &
1034        ment, elij, delt, plcl, mp, qp, up, vp, trap, wt, water, evap, b, &
1035        wdtraina, wdtrainm &
1036#ifdef ISO
1037     &               ,xt,xtclw,xtelij &
1038     &               ,xtp,xtwater,xtevap &
1039#ifdef DIAGISO
1040     &               ,wdtrain,xtwdtrain,taux_cond_conv &
1041#endif       
1042#endif
1043     &               )
1044      ! RomP <<<
1045
1046#ifdef ISO
1047!       write(*,*) 'klev=',klev
1048#ifdef ISOVERIF
1049       write(*,*) 'cv_driver 930: après cv3_unsat'
1050       do k = 1, klev
1051        do i = 1, ncum
1052         if (use_iso(iso_eau)) then
1053            call iso_verif_egalite_choix(xt(iso_eau,i,k),q(i,k), &
1054     &            'cv_driver 778',errmax,errmaxrel)
1055            call iso_verif_egalite_choix(xtp(iso_eau,i,k), &
1056     &            qp(i,k),'cv_driver 780',errmax,errmaxrel)
1057            call iso_verif_egalite_choix( &
1058     &             xtwater(iso_eau,i,k),water(i,k), &
1059     &            'cv_driver 782',errmax,errmaxrel)
1060            call iso_verif_egalite_choix(xtevap(iso_eau,i,k), &
1061     &            evap(i,k),'cv_driver 784',errmax,errmaxrel)
1062         endif !  if (use_iso(iso_eau)) then
1063#ifdef ISOTRAC
1064           call iso_verif_traceur(xt(1,i,k),'cv_driver 939')
1065           call iso_verif_traceur(xtp(1,i,k),'cv_driver 940')
1066           call iso_verif_traceur(xtwater(1,i,k),'cv_driver 941')
1067           call iso_verif_traceur_justmass(xtevap(1,i,k), &
1068     &           'cv_driver 942')
1069           if (bidouille_anti_divergence) then
1070             call iso_verif_traceur_pbidouille( &
1071     &           xtp(1,i,k),'cv_driver 956')
1072             call iso_verif_traceur_pbidouille( &
1073     &           xtwater(1,i,k),'cv_driver 958')
1074           endif
1075#endif             
1076        enddo !do i = 1, ncum
1077       enddo !do k = 1, klev     
1078#endif
1079#endif
1080
1081    END IF
1082
1083    IF (iflag_con==4) THEN
1084      CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph, h, lv, &
1085        ep, sigp, clw, m, ment, elij, iflag, mp, qp, up, vp, wt, water, evap)
1086    END IF
1087
1088#ifdef ISO
1089#ifdef ISOTRAC
1090      if (option_traceurs.eq.6) then
1091          ! on colorie les ddfts en rouge, le reste est en blanc.
1092          do k = 1, klev
1093            do i = 1, ncum
1094               do iiso=1,niso
1095                  ixt_ddft=index_trac(izone_ddft,iiso)
1096                  ixt_poubelle=index_trac(izone_poubelle,iiso)
1097                  xtp(ixt_ddft,i,k)=xtp(ixt_ddft,i,k) &
1098     &                    +xtp(ixt_poubelle,i,k)
1099                  xtp(ixt_poubelle,i,k)=0.0
1100               enddo !do iiso=1,niso
1101#ifdef ISOVERIF
1102               call iso_verif_traceur(xtp(1,i,k),'cv_driver 990')
1103#endif               
1104            enddo !do i = 1, ncum
1105          enddo !do k = 1, klev
1106      else if (option_traceurs.eq.19) then
1107          ! on colorie les ddfts, mais on sauve la revap
1108          do k = 1, klev
1109            do i = 1, ncum
1110               do izone=1,ntraceurs_zone
1111                 if (izone.eq.izone_ddft) then
1112                   do iiso=1,niso
1113                     ixt_ddft=index_trac(izone,iiso)
1114                     ixt_revap=index_trac(izone_revap,iiso)
1115                     xtp(ixt_ddft,i,k)=xtp(iiso,i,k)-xtp(ixt_revap,i,k)
1116                   enddo !do iiso=1,niso
1117                 elseif (izone.eq.izone_ddft) then
1118                    ! rien à faire
1119                 else !if (izone.eq.izone_ddft) then
1120                   do iiso=1,niso
1121                     ixt=index_trac(izone,iiso)
1122                     xtp(ixt,i,k)=0.0
1123                   enddo !do iiso=1,niso
1124                 endif !if (izone.eq.izone_ddft) then
1125               enddo !do izone=1,ntraceurs_zone
1126#ifdef ISOVERIF
1127               call iso_verif_traceur(xtp(1,i,k),'cv_driver 1059')
1128#endif               
1129            enddo !do i = 1, ncum
1130          enddo !do k = 1, klev
1131      endif !if (option_traceurs.eq.6) then
1132#endif
1133#endif   
1134    ! -------------------------------------------------------------------
1135    ! --- YIELD
1136    ! (tendencies, precipitation, variables of interface with other
1137    ! processes, etc)
1138    ! -------------------------------------------------------------------
1139
1140    IF (iflag_con==30) THEN
1141      CALL cv30_yield(nloc, ncum, nd, nd, ntra & ! na->nd
1142        , icb, inb, delt, t, q, u, v, tra, gz, p, ph, h, hp, lv, cpn, th, ep, &
1143        clw, m, tp, mp, qp, up, vp, trap, wt, water, evap, b, ment, qent, &
1144        uent, vent, nent, elij, traent, sig, tv, tvp, iflag, precip, vprecip, &
1145        ft, fq, fu, fv, ftra, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc, &
1146        wd &
1147#ifdef ISO
1148     &                     ,xt,xtclw,xtp,xtwater,xtevap &
1149     &                     ,xtent,xtelij,xtprecip,fxt,xtVPrecip &
1150#ifdef DIAGISO
1151     &         ,fq_detrainement,fq_ddft,fq_fluxmasse,fq_evapprecip  &
1152     &         ,fxt_detrainement,fxt_ddft,fxt_fluxmasse,fxt_evapprecip &
1153     &         ,Amp_diag     &
1154     &         ,f_detrainement,q_detrainement,xt_detrainement  &
1155#endif       
1156#endif
1157     &      )
1158
1159#ifdef ISOVERIF
1160      DO k = 1, nd
1161       DO i = 1, ncum     
1162        if (use_iso(iso_eau)) then
1163             call iso_verif_egalite_choix(fxt(iso_eau,i,k), &
1164     &           fq(i,k),'cv_driver 1086',errmax,errmaxrel)
1165        endif !if (use_iso(iso_HDO)) then
1166        if (use_iso(iso_HDO)) then
1167          if (q(i,k).gt.ridicule) then
1168            call iso_verif_aberrant( &
1169     &          (xt(iso_HDO,i,k)+delt*fxt(iso_HDO,i,k)) &
1170     &          /(q(i,k)+delt*fq(i,k)),'cv_driver 855')
1171          endif
1172         endif
1173#ifdef ISOTRAC
1174           call iso_verif_traceur(xt(1,i,k),'cv_driver 1008')
1175#endif       
1176#ifdef DIAGISO
1177        call iso_verif_positif(tg_save(i,k)-50.0,'cv_driver 1024')
1178#endif             
1179        enddo
1180       enddo
1181#endif
1182
1183    END IF
1184
1185    IF (iflag_con==4) THEN
1186      CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt, t, q, u, v, gz, p, &
1187        ph, h, hp, lv, cpn, ep, clw, frac, m, mp, qp, up, vp, wt, water, &
1188        evap, ment, qent, uent, vent, nent, elij, tv, tvp, iflag, wd, qprime, &
1189        tprime, precip, cbmf, ft, fq, fu, fv, ma, qcondc)
1190    END IF
1191
1192    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1193    ! --- passive tracers
1194    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1195
1196    IF (iflag_con==30) THEN
1197      ! RomP >>>
1198      CALL cv30_tracer(nloc, len, ncum, nd, nd, ment, sij, da, phi, phi2, &
1199        d1a, dam, ep, vprecip, elij, clw, epmlmmm, eplamm, icb, inb)
1200      ! RomP <<<
1201    END IF
1202
1203    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1204    ! --- UNCOMPRESS THE FIELDS
1205    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1206    ! set iflag1 =42 for non convective points
1207    DO i = 1, len
1208      iflag1(i) = 42
1209    END DO
1210
1211    IF (iflag_con==30) THEN
1212      CALL cv30_uncompress(nloc, len, ncum, nd, ntra, idcum, iflag, precip, &
1213        vprecip, evap, ep, sig, w0 & !RomP
1214        , ft, fq, fu, fv, ftra, inb, ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
1215        da, phi, mp, phi2, d1a, dam, sij & !RomP
1216        , elij, clw, epmlmmm, eplamm & !RomP
1217        , wdtraina, wdtrainm,epmax_diag &     !RomP
1218        , iflag1, precip1, vprecip1, evap1, ep1, sig1, w01 & !RomP
1219        , ft1, fq1, fu1, fv1, ftra1, inb1, ma1, upwd1, dnwd1, dnwd01, &
1220        qcondc1, wd1, cape1, da1, phi1, mp1, phi21, d1a1, dam1, sij1 & !RomP
1221        , elij1, clw1, epmlmmm1, eplamm1 & !RomP
1222        , wdtraina1, wdtrainm1,epmax_diag1 & !RomP
1223#ifdef ISO
1224     &          ,xtprecip,fxt,xtVPrecip   & 
1225     &         ,xtprecip1,fxt1,xtVPrecip1 &
1226#ifdef DIAGISO
1227     &         , water,xtwater,qp,xtp,xtevap &
1228     &         , xtclw &
1229     &         , wdtrain,xtwdtrain,taux_cond_conv &
1230     &         , fq_detrainement,fq_ddft,fq_fluxmasse,fq_evapprecip &
1231     &         , fxt_detrainement,fxt_ddft,fxt_fluxmasse &
1232     &         , fxt_evapprecip,m,Amp_diag,tg_save,epmax_diag &
1233     &         , f_detrainement,q_detrainement,xt_detrainement &
1234!     :         , ment,qent,xtent
1235     &         , water1,xtwater1,qp1,xtp1,xtevap1 &
1236     &         , xtclw1 &
1237     &         , wdtrain1,xtwdtrain1,taux_cond_conv1 &
1238     &         , fq_detrainement1,fq_ddft1,fq_fluxmasse1,fq_evapprecip1  &
1239     &         , fxt_detrainement1,fxt_ddft1,fxt_fluxmasse1 &
1240     &         , fxt_evapprecip1,Mi1,Amp_diag1,tg_save1,epmax_diag1 &
1241     &         , f_detrainement1,q_detrainement1,xt_detrainement1 &
1242#endif       
1243#endif
1244     &          )
1245
1246#ifdef ISO       
1247#ifdef ISOVERIF
1248      DO k = 1, nd
1249       DO i = 1, len       
1250        if (use_iso(iso_eau)) then
1251             call iso_verif_egalite_choix(fxt1(iso_eau,i,k), &
1252     &          fq1(i,k),'cv_driver 1170',errmax,errmaxrel)
1253        endif !if (use_iso(iso_HDO)) then
1254        if (use_iso(iso_HDO)) then
1255          if (q1(i,k).gt.ridicule) then
1256            call iso_verif_aberrant( &
1257     &          (xt1(iso_HDO,i,k)+delt*fxt1(iso_HDO,i,k)) &
1258     &          /(q1(i,k)+delt*fq1(i,k)),'cv_driver 2554')
1259          endif
1260         endif         
1261#ifdef DIAGISO
1262         call iso_verif_positif(tg_save1(i,k)-50.0,'cv_driver 1221')
1263#endif         
1264        enddo !DO i = 1, len
1265       enddo !DO k = 1, nd
1266#ifdef ISOTRAC       
1267           call iso_verif_trac_masse_vect(fxt1,len,nd, &
1268     &           'cv_driver 1200',errmax,errmaxrel)   
1269#endif             
1270#endif
1271#endif
1272
1273    END IF
1274
1275    IF (iflag_con==4) THEN
1276      CALL cv_uncompress(nloc, len, ncum, nd, idcum, iflag, precip, cbmf, ft, &
1277        fq, fu, fv, ma, qcondc, iflag1, precip1, cbmf1, ft1, fq1, fu1, fv1, &
1278        ma1, qcondc1)
1279    END IF
1280
1281  END IF ! ncum>0
1282
1283  ! print *, 'fin cv_driver ->'      !jyg
1284  RETURN
1285END SUBROUTINE cv_driver
1286
1287! ==================================================================
1288SUBROUTINE cv_flag(iflag_ice_thermo)
1289  IMPLICIT NONE
1290
1291  ! Argument : iflag_ice_thermo : ice thermodynamics is taken into account if
1292  ! iflag_ice_thermo >=1
1293  INTEGER iflag_ice_thermo
1294
1295  include "cvflag.h"
1296
1297  ! -- si .TRUE., on rend la gravite plus explicite et eventuellement
1298  ! differente de 10.0 dans convect3:
1299  cvflag_grav = .TRUE.
1300  cvflag_ice = iflag_ice_thermo >= 1
1301
1302  RETURN
1303END SUBROUTINE cv_flag
1304
1305! ==================================================================
1306SUBROUTINE cv_thermo(iflag_con)
1307  IMPLICIT NONE
1308
1309  ! -------------------------------------------------------------
1310  ! Set thermodynamical constants for convectL
1311  ! -------------------------------------------------------------
1312
1313  include "YOMCST.h"
1314  include "cvthermo.h"
1315
1316  INTEGER iflag_con
1317
1318
1319  ! original set from convect:
1320  IF (iflag_con==4) THEN
1321    cpd = 1005.7
1322    cpv = 1870.0
1323    cl = 4190.0
1324    rrv = 461.5
1325    rrd = 287.04
1326    lv0 = 2.501E6
1327    g = 9.8
1328    t0 = 273.15
1329    grav = g
1330  ELSE
1331
1332    ! constants consistent with LMDZ:
1333    cpd = rcpd
1334    cpv = rcpv
1335    cl = rcw
1336    ci = rcs
1337    rrv = rv
1338    rrd = rd
1339    lv0 = rlvtt
1340    lf0 = rlstt - rlvtt
1341    g = rg ! not used in convect3
1342    ! ori      t0  = RTT
1343    t0 = 273.15 ! convect3 (RTT=273.16)
1344    ! maf       grav= 10.    ! implicitely or explicitely used in convect3
1345    grav = g ! implicitely or explicitely used in convect3
1346  END IF
1347
1348  rowl = 1000.0 !(a quelle variable de YOMCST cela correspond-il?)
1349
1350  clmcpv = cl - cpv
1351  clmcpd = cl - cpd
1352  clmci = cl - ci
1353  cpdmcp = cpd - cpv
1354  cpvmcpd = cpv - cpd
1355  cpvmcl = cl - cpv ! for convect3
1356  eps = rrd/rrv
1357  epsi = 1.0/eps
1358  epsim1 = epsi - 1.0
1359  ! ginv=1.0/g
1360  ginv = 1.0/grav
1361  hrd = 0.5*rrd
1362
1363  RETURN
1364END SUBROUTINE cv_thermo
Note: See TracBrowser for help on using the repository browser.