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

Last change on this file since 5209 was 5160, checked in by abarral, 7 weeks ago

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