source: LMDZ6/branches/Amaury_dev/libf/phylmdiso/cv_driver.F90 @ 5141

Last change on this file since 5141 was 5141, checked in by abarral, 3 months ago

Put cvthermo.h, cv30param.h, cv3param.h into modules

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