source: LMDZ6/trunk/libf/phylmdiso/cv_driver.F90 @ 4738

Last change on this file since 4738 was 4491, checked in by crisi, 20 months ago

Bug corrections in LMDZiso, especially for water tagging

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