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

Last change on this file since 5435 was 5285, checked in by abarral, 8 weeks ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

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