module cv3a_driver_mod implicit none contains ! ------------------------------------------------------------------- ! Prolog by Kerry Emanuel. ! ------------------------------------------------------------------- ! --- ARGUMENTS ! ------------------------------------------------------------------- ! --- On input: ! ! t: Array of absolute temperature (K) of dimension ND, with first ! index corresponding to lowest model level. Note that this array ! will be altered by the subroutine if dry convective adjustment ! occurs and if IPBL is not equal to 0. ! ! q: Array of specific humidity (gm/gm) of dimension ND, with first ! index corresponding to lowest model level. Must be defined ! at same grid levels as T. Note that this array will be altered ! if dry convective adjustment occurs and if IPBL is not equal to 0. ! ! qs: Array of saturation specific humidity of dimension ND, with first ! index corresponding to lowest model level. Must be defined ! at same grid levels as T. Note that this array will be altered ! if dry convective adjustment occurs and if IPBL is not equal to 0. ! ! t_wake: Array of absolute temperature (K), seen by unsaturated draughts, ! of dimension ND, with first index corresponding to lowest model level. ! ! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts, ! of dimension ND, with first index corresponding to lowest model level. ! Must be defined at same grid levels as T. ! ! qs_wake: Array of saturation specific humidity, seen by unsaturated draughts, ! of dimension ND, with first index corresponding to lowest model level. ! Must be defined at same grid levels as T. ! ! s_wake: Array of fractionnal area occupied by the wakes. ! ! u: Array of zonal wind velocity (m/s) of dimension ND, witth first ! index corresponding with the lowest model level. Defined at ! same levels as T. Note that this array will be altered if ! dry convective adjustment occurs and if IPBL is not equal to 0. ! ! v: Same as u but for meridional velocity. ! ! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA), ! where NTRA is the number of different tracers. If no ! convective tracer transport is needed, define a dummy ! input array of dimension (ND,1). Tracers are defined at ! same vertical levels as T. Note that this array will be altered ! if dry convective adjustment occurs and if IPBL is not equal to 0. ! ! p: Array of pressure (mb) of dimension ND, with first ! index corresponding to lowest model level. Must be defined ! at same grid levels as T. ! ! ph: Array of pressure (mb) of dimension ND+1, with first index ! corresponding to lowest level. These pressures are defined at ! levels intermediate between those of P, T, Q and QS. The first ! value of PH should be greater than (i.e. at a lower level than) ! the first value of the array P. ! ! ALE: Available lifting Energy ! ! ALP: Available lifting Power ! ! nl: The maximum number of levels to which convection can penetrate, plus 1. ! NL MUST be less than or equal to ND-1. ! ! delt: The model time step (sec) between calls to CONVECT ! ! ---------------------------------------------------------------------------- ! --- On Output: ! ! iflag: An output integer whose value denotes the following: ! VALUE INTERPRETATION ! ----- -------------- ! 0 Moist convection occurs. ! 1 Moist convection occurs, but a CFL condition ! on the subsidence warming is violated. This ! does not cause the scheme to terminate. ! 2 Moist convection, but no precip because ep(inb) lt 0.0001 ! 3 No moist convection because new cbmf is 0 and old cbmf is 0. ! 4 No moist convection; atmosphere is not ! unstable ! 6 No moist convection because ihmin le minorig. ! 7 No moist convection because unreasonable ! parcel level temperature or specific humidity. ! 8 No moist convection: lifted condensation ! level is above the 200 mb level. ! 9 No moist convection: cloud base is higher ! then the level NL-1. ! 10 No moist convection: cloud top is too warm. ! ! ! ft: Array of temperature tendency (K/s) of dimension ND, defined at same ! grid levels as T, Q, QS and P. ! ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND, ! defined at same grid levels as T, Q, QS and P. ! ! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND, ! defined at same grid levels as T. ! ! fv: Same as FU, but for forcing of meridional velocity. ! ! ftra: Array of forcing of tracer content, in tracer mixing ratio per ! second, defined at same levels as T. Dimensioned (ND,NTRA). ! ! precip: Scalar convective precipitation rate (mm/day). ! ! wd: A convective downdraft velocity scale. For use in surface ! flux parameterizations. See convect.ps file for details. ! ! tprime: A convective downdraft temperature perturbation scale (K). ! For use in surface flux parameterizations. See convect.ps ! file for details. ! ! qprime: A convective downdraft specific humidity ! perturbation scale (gm/gm). ! For use in surface flux parameterizations. See convect.ps ! file for details. ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT ! ITS NEXT CALL. That is, the value of CBMF must be "remembered" ! by the calling program between calls to CONVECT. ! ! det: Array of detrainment mass flux of dimension ND. ! ------------------------------------------------------------------- ! written by : Sandrine Bony-Lena , 17/05/2003, 11.19.41 ! S. Bony, Mar 2002: ! * Several modules corresponding to different physical processes ! * Several versions of convect may be used: ! - iflag_con=3: version lmd (previously named convect3) ! - iflag_con=4: version 4.3b (vect. version, previously convect1/2) ! + tard: - iflag_con=5: version lmd with ice (previously named convectg) ! S. Bony, Oct 2002: ! * Vectorization of convect3 (ie version lmd) SUBROUTINE cv3a_driver(len, nd, ndp1, ntra, nloc, k_upper, & iflag_con, iflag_mix, iflag_ice_thermo, iflag_clos, ok_conserv_q, & delt, comp_threshold, & t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, & u1, v1, & p1, ph1, & Ale1, Alp1, omega1, & sig1feed1, sig2feed1, wght1, & iflag1, ft1, fq1, fu1, fv1, ftra1, & precip1, kbas1, ktop1, & cbmf1, plcl1, plfc1, wbeff1, & sig1, w01, & !input/output ptop21, sigd1, & ma1, mip1, Vprecip1, Vprecipi1, upwd1, dnwd1, dnwd01, & qcondc1, wd1, & cape1, cin1, tvp1, & ftd1, fqd1, & Plim11, Plim21, asupmax1, supmax01, asupmaxmin1, & da1, phi1, mp1, phi21, d1a1, dam1, sigij1, wghti1, & qta1, clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, & wdtrainA1, wdtrainS1, wdtrainM1, qtc1, sigt1, tau_cld_cv, & coefw_cld_cv, & epmax_diag1) USE print_control_mod, ONLY: prt_level, lunout USE cv3a_compress_mod USE profiling_physic_mod IMPLICIT NONE ! Input INTEGER, INTENT(IN) :: len ! first (i) dimension INTEGER, INTENT(IN) :: nd ! vertical (k) dimension INTEGER, INTENT(IN) :: ndp1 ! nd + 1 INTEGER, INTENT(IN) :: ntra ! number of tracors INTEGER, INTENT(IN) :: k_upper ! upmost level for vertical loops INTEGER, INTENT(IN) :: iflag_con ! version of convect (3) INTEGER, INTENT(IN) :: iflag_mix ! version of mixing (0/1/2) INTEGER, INTENT(IN) :: iflag_ice_thermo ! accounting for ice thermodynamics (0/1) INTEGER, INTENT(IN) :: iflag_clos ! version of closure (0/1) LOGICAL, INTENT(IN) :: ok_conserv_q ! when true corrections for water conservation are swtiched on REAL, INTENT(IN) :: tau_cld_cv ! characteristic time of dissipation of mixing fluxes REAL, INTENT(IN) :: coefw_cld_cv ! coefficient for updraft velocity in convection REAL, INTENT(IN) :: delt ! time step REAL, INTENT (IN) :: comp_threshold REAL, DIMENSION(len, nd), INTENT(IN) :: t1 ! temperature (sat draught envt) REAL, DIMENSION(len, nd), INTENT(IN) :: q1 ! specific hum (sat draught envt) REAL, DIMENSION(len, nd), INTENT(IN) :: qs1 ! sat specific hum (sat draught envt) REAL, DIMENSION(len, nd), INTENT(IN) :: t1_wake ! temperature (unsat draught envt) REAL, DIMENSION(len, nd), INTENT(IN) :: q1_wake ! specific hum(unsat draught envt) REAL, DIMENSION(len, nd), INTENT(IN) :: qs1_wake ! sat specific hum(unsat draughts envt) REAL, DIMENSION(len), INTENT(IN) :: s1_wake ! fractionnal area covered by wakes REAL, DIMENSION(len, nd), INTENT(IN) :: u1 ! u-wind REAL, DIMENSION(len, nd), INTENT(IN) :: v1 ! v-wind REAL, DIMENSION(len, nd), INTENT(IN) :: p1 ! full level pressure REAL, DIMENSION(len, ndp1), INTENT(IN) :: ph1 ! half level pressure REAL, DIMENSION(len), INTENT(IN) :: Ale1 ! Available lifting Energy REAL, DIMENSION(len), INTENT(IN) :: Alp1 ! Available lifting Power REAL, DIMENSION(len, nd), INTENT(IN) :: omega1 REAL, INTENT(IN) :: sig1feed1 ! sigma coord/pressure at lower bound of feeding layer REAL, INTENT(IN) :: sig2feed1 ! sigma coord/pressure at upper bound of feeding layer REAL, DIMENSION(nd), INTENT(IN) :: wght1 ! weight density determining the feeding mixture ! Input/Output REAL, DIMENSION(len, nd), INTENT(INOUT) :: sig1 ! section adiabatic updraft REAL, DIMENSION(len, nd), INTENT(INOUT) :: w01 ! vertical velocity within adiab updraft ! Output INTEGER, DIMENSION(len), INTENT(OUT) :: iflag1 ! flag for Emanuel conditions REAL, DIMENSION(len, nd), INTENT(OUT) :: ft1 ! temp tend REAL, DIMENSION(len, nd), INTENT(OUT) :: fq1 ! spec hum tend REAL, DIMENSION(len, nd), INTENT(OUT) :: fu1 ! u-wind tend REAL, DIMENSION(len, nd), INTENT(OUT) :: fv1 ! v-wind tend REAL, DIMENSION(len, nd, ntra), INTENT(OUT) :: ftra1 ! tracor tend REAL, DIMENSION(len), INTENT(OUT) :: precip1 ! precipitation INTEGER, DIMENSION(len), INTENT(OUT) :: kbas1 ! cloud base level INTEGER, DIMENSION(len), INTENT(OUT) :: ktop1 ! cloud top level REAL, DIMENSION(len), INTENT(OUT) :: cbmf1 ! cloud base mass flux REAL, DIMENSION(len), INTENT(OUT) :: plcl1 REAL, DIMENSION(len), INTENT(OUT) :: plfc1 REAL, DIMENSION(len), INTENT(OUT) :: wbeff1 REAL, DIMENSION(len), INTENT(OUT) :: ptop21 ! top of entraining zone REAL, DIMENSION(len), INTENT(OUT) :: sigd1 REAL, DIMENSION(len, nd), INTENT(OUT) :: ma1 ! mass flux adiabatic updraft (staggered grid) REAL, DIMENSION(len, nd), INTENT(OUT) :: mip1 ! mass flux shed by the adiabatic updraft (extensive) REAL, DIMENSION(len, ndp1), INTENT(OUT) :: vprecip1 ! vertical profile of total precipitation (staggered grid) REAL, DIMENSION(len, ndp1), INTENT(OUT) :: vprecipi1 ! vertical profile of ice precipitation (staggered grid) REAL, DIMENSION(len, nd), INTENT(OUT) :: upwd1 ! total upward mass flux (adiab+mixed) (staggered grid) REAL, DIMENSION(len, nd), INTENT(OUT) :: dnwd1 ! saturated downward mass flux (mixed) (staggered grid) REAL, DIMENSION(len, nd), INTENT(OUT) :: dnwd01 ! unsaturated downward mass flux (staggered grid) REAL, DIMENSION(len, nd), INTENT(OUT) :: qcondc1 ! in-cld mixing ratio of condensed water (intensive) REAL, DIMENSION(len), INTENT(OUT) :: wd1 ! downdraft velocity scale for sfc fluxes REAL, DIMENSION(len), INTENT(OUT) :: cape1 REAL, DIMENSION(len), INTENT(OUT) :: cin1 REAL, DIMENSION(len, nd), INTENT(OUT) :: tvp1 ! adiab lifted parcell virt temp REAL, DIMENSION(len, nd), INTENT(OUT) :: ftd1 ! temperature tendency due to precipitations (K/s) REAL, DIMENSION(len, nd), INTENT(OUT) :: fqd1 ! specific humidity tendencies due to precipitations ((gm/gm)/s) REAL, DIMENSION(len), INTENT(OUT) :: Plim11 REAL, DIMENSION(len), INTENT(OUT) :: Plim21 REAL, DIMENSION(len, nd), INTENT(OUT) :: asupmax1 ! Highest mixing fraction of mixed updraughts REAL, DIMENSION(len), INTENT(OUT) :: supmax01 REAL, DIMENSION(len), INTENT(OUT) :: asupmaxmin1 REAL, DIMENSION(len, nd), INTENT(OUT) :: qtc1 ! in cloud water content / specific humidity in convection (intensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: sigt1 ! surface fraction in adiabatic updrafts / fract. cloud area (intensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: wdtrainA1 ! precipitation ejected from adiabatic draught REAL, DIMENSION(len, nd), INTENT(OUT) :: wdtrainS1 ! precipitation detrained from shedding of adiabatic draught REAL, DIMENSION(len, nd), INTENT(OUT) :: wdtrainM1 ! precipitation detrained from mixed draughts REAL, DIMENSION(len, nd), INTENT(OUT) :: mp1 ! unsat. mass flux (staggered grid) REAL, DIMENSION(len, nd), INTENT(OUT) :: da1 ! detrained mass flux of adiab. asc. air (extensive) REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: phi1 ! mass flux of envt. air in mixed draughts (extensive) REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: epmlmMm1 ! (extensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: eplaMm1 ! (extensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: evap1 ! evaporation rate in precip. downdraft. (intensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: ep1 REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: sigij1 ! mass fraction of env. air in mixed draughts (intensive) REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: elij1! cond. water per unit mass of mixed draughts (intensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: qta1 ! total water per unit mass of the adiab. asc. (intensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: clw1 ! condensed water content of the adiabatic updraught / cond. water per unit mass of the adiab. asc. (intensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: wghti1 ! final weight of the feeding layers (extensive) REAL, DIMENSION(len, nd, nd), INTENT(OUT) :: phi21 ! (extensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: d1a1 ! (extensive) REAL, DIMENSION(len, nd), INTENT(OUT) :: dam1 ! (extensive) REAL, DIMENSION(len), INTENT(OUT) :: epmax_diag1 ! Local (non compressed) arrays INTEGER nk1(len) INTEGER icbs1(len) LOGICAL, SAVE :: debut = .TRUE. !$omp THREADPRIVATE(debut) REAL tnk1(len) REAL qnk1(len) REAL gznk1(len) REAL unk1(len) REAL vnk1(len) REAL hnk1(len) REAL pbase1(len) REAL buoybase1(len) REAL lf1(len, nd), lf1_wake(len, nd) REAL lv1(len, nd), lv1_wake(len, nd) REAL cpn1(len, nd), cpn1_wake(len, nd) REAL tv1(len, nd), tv1_wake(len, nd) REAL gz1(len, nd), gz1_wake(len, nd) REAL h1(len, nd), h1_wake(len, nd) REAL tp1(len, nd) REAL th1(len, nd), th1_wake(len, nd) CHARACTER(LEN=20) :: modname = 'cva_driver' CHARACTER(LEN=80) :: abort_message INTEGER :: coef_convective(len) INTEGER :: k REAL, parameter :: Cin_noconv = -100000., Cape_noconv = -1 INTEGER, SAVE :: igout = 1 !$omp THREADPRIVATE(igout) INTEGER nloc ! Allocation size for compressed arrays logical :: compress call enter_profile("cv3a_driver") if (iflag_con /= 3) call abort_physic("cv3a_driver", "iflag_con must be 3", 1) ! ------------------------------------------------------------------- ! --- SET CONSTANTS AND PARAMETERS ! ------------------------------------------------------------------- ! -- set simulation flags: ! (common cvflag CALL cv_flag(iflag_ice_thermo) ! -- set thermodynamical constants: ! (common cvthermo) CALL cv_thermo(iflag_con) ! -- set convect parameters ! includes microphysical parameters and parameters that ! control the rate of approach to quasi-equilibrium) ! (common cvparam) CALL cv3_param(nd, k_upper, delt) CALL cv3_incrcount(len, nd, delt, sig1) ! -------------------------------------------------------------------- ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY ! -------------------------------------------------------------------- IF (debut) PRINT *, 'Emanuel version 3 nouvelle' block REAL thnk1(len) REAL qsnk1(len) REAL cpnk1(len) REAL hm1(len, nd) REAL bid(len, nd) ! dummy array REAL p1feed1(len) ! pressure at lower bound of feeding layer REAL p2feed1(len) ! pressure at upper bound of feeding layer integer :: i, icbmax call driver_log('cv3_prelim') CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, & lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1) call driver_log('cv3_prelim') CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, & lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, & h1_wake, bid, th1_wake) ! -------------------------------------------------------------------- ! --- CONVECTIVE FEED ! -------------------------------------------------------------------- ! compute feeding layer potential temperature and mixing ratio : ! get bounds of feeding layer ! test niveaux couche alimentation KE IF (sig1feed1 == sig2feed1) THEN WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed' WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def' abort_message = '' CALL abort_physic(modname, abort_message, 1) END IF DO i = 1, len p1feed1(i) = sig1feed1*ph1(i, 1) p2feed1(i) = sig2feed1*ph1(i, 1) END DO call driver_log('cv3_feed') ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed iflag1(:) = 0 plcl1(:) = 0. wghti1(:, :) = 0. CALL cv3_feed(len, nd, ok_conserv_q, & t1, q1, u1, v1, p1, ph1, h1, gz1, & p1feed1, p2feed1, wght1, & wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, & cpnk1, hnk1, nk1, kbas1, icbmax, iflag1, gznk1, plcl1) ! -------------------------------------------------------------------- ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part ! (up through ICB for convect4, up through ICB+1 for convect3) ! Calculates the lifted parcel virtual temperature at nk, the ! actual temperature, and the adiabatic liquid water content. ! -------------------------------------------------------------------- call driver_log('cv3_undilute1') ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed tvp1(:, :) = 0. clw1(:, :) = 0. CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, kbas1, tnk1, qnk1, & ! nd->na gznk1, tp1, tvp1, clw1, icbs1) ! ------------------------------------------------------------------- ! --- TRIGGERING ! ------------------------------------------------------------------- call driver_log('cv3_trigger') CALL cv3_trigger(len, nd, kbas1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na pbase1, buoybase1, iflag1, sig1, w01) end block ! ===================================================================== ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY ! ===================================================================== ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ ! --- COMPRESS THE FIELDS ! (-> vectorization over convective gridpoints) ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ ! Compression has 3 possible modes: ! 1) Compress = true : compress_mode = COMPRESS_MODE_COMPRESS ! 1.1) Copy convective cells in contiguous array !compress_mode = COMPRESS_MODE_COPY ! 1.2) Copy all cells and also compute on non-convective cells ! 2) Compress = false : Don't copy and use original arrays, compute on non-convective cells ! Never compress when comp_threshold = 0, always compress when comp_threshold = 1 ! comp_threshold = 0 (default) offers best performance when using many CPUs compress = (comp_threshold == 1) .or. ( get_compress_size(len, (iflag1(:) == 0)) < len*comp_threshold ) if (compress) then nloc = get_compress_size(len, (iflag1(:) == 0)) block ! (local) compressed fields: INTEGER iflag(nloc), nk(nloc), icb(nloc) INTEGER nent(nloc, nd) INTEGER icbs(nloc) INTEGER inb(nloc) REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc) REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd) REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd) REAL s_wake(nloc) REAL u(nloc, nd), v(nloc, nd) REAL gz(nloc, nd), h(nloc, nd) REAL h_wake(nloc, nd) REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd) REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd) REAL p(nloc, nd), ph(nloc, nd + 1), tv(nloc, nd), tp(nloc, nd) REAL tv_wake(nloc, nd) REAL clw(nloc, nd) REAL qta(nloc, nd), qpreca(nloc, nd) REAL pbase(nloc), buoybase(nloc), th(nloc, nd) REAL th_wake(nloc, nd) REAL tvp(nloc, nd) REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc) REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd) REAL buoy(nloc, nd) REAL cape(nloc) REAL cin(nloc) REAL m(nloc, nd) REAL mm(nloc, nd) REAL ment(nloc, nd, nd), sigij(nloc, nd, nd) REAL qent(nloc, nd, nd) REAL hent(nloc, nd, nd) REAL uent(nloc, nd, nd), vent(nloc, nd, nd) REAL ments(nloc, nd, nd), qents(nloc, nd, nd) REAL elij(nloc, nd, nd) REAL supmax(nloc, nd) REAL Ale(nloc), Alp(nloc), coef_clos(nloc) REAL omega(nloc, nd) REAL sigd(nloc) REAL, DIMENSION(nloc, nd) :: mp, qp, up, vp REAL, DIMENSION(nloc, nd) :: wt, water, evap REAL, DIMENSION(nloc, nd) :: ice, fondue, b REAL, DIMENSION(nloc, nd) :: frac_a, frac_s, faci REAL ft(nloc, nd), fq(nloc, nd) REAL ftd(nloc, nd), fqd(nloc, nd) REAL fu(nloc, nd), fv(nloc, nd) REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd) REAL ma(nloc, nd), mip(nloc, nd) REAL precip(nloc) REAL vprecip(nloc, nd + 1) REAL vprecipi(nloc, nd + 1) REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra) REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra) REAL qcondc(nloc, nd) REAL wd(nloc) REAL Plim1(nloc), plim2(nloc) REAL asupmax(nloc, nd) REAL supmax0(nloc) REAL asupmaxmin(nloc) REAL tnk(nloc), qnk(nloc), gznk(nloc) REAL wghti(nloc, nd) REAL hnk(nloc), unk(nloc), vnk(nloc) REAL qtc(nloc, nd) REAL sigt(nloc, nd) REAL wdtrainA(nloc, nd), wdtrainS(nloc, nd), wdtrainM(nloc, nd) REAL da(nloc, nd), phi(nloc, nd, nd) REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd) REAL phi2(nloc, nd, nd) REAL d1a(nloc, nd), dam(nloc, nd) REAL epmax_diag(nloc) LOGICAL ok_inhib ! True => possible inhibition of convection by dryness integer :: ncum ! Number of convective cells type(compress_data_t) :: compress_data type(array_list) :: cv3a_compress_list, cv3a_uncompress_list logical, save :: timers_first = .true. !$omp threadprivate(timers_first) call enter_profile("cv3a_compress") call driver_log('cv3a_compress') call add_array_i1(cv3a_compress_list, iflag1, iflag) call add_array_i1(cv3a_compress_list, nk1, nk) call add_array_i1(cv3a_compress_list, kbas1, icb) call add_array_i1(cv3a_compress_list, icbs1, icbs) call add_array_r1(cv3a_compress_list, plcl1, plcl) call add_array_r1(cv3a_compress_list, tnk1, tnk) call add_array_r1(cv3a_compress_list, qnk1, qnk) call add_array_r1(cv3a_compress_list, gznk1, gznk) call add_array_r1(cv3a_compress_list, hnk1, hnk) call add_array_r1(cv3a_compress_list, unk1, unk) call add_array_r1(cv3a_compress_list, vnk1, vnk) call add_array_r2(cv3a_compress_list, wghti1, wghti) call add_array_r1(cv3a_compress_list, pbase1, pbase) call add_array_r1(cv3a_compress_list, buoybase1, buoybase) call add_array_r2(cv3a_compress_list, th1, th) call add_array_r2(cv3a_compress_list, t1, t) call add_array_r2(cv3a_compress_list, q1, q) call add_array_r2(cv3a_compress_list, qs1, qs) call add_array_r2(cv3a_compress_list, t1_wake, t_wake) call add_array_r2(cv3a_compress_list, q1_wake, q_wake) call add_array_r2(cv3a_compress_list, qs1_wake, qs_wake) call add_array_r1(cv3a_compress_list, s1_wake, s_wake) call add_array_r2(cv3a_compress_list, u1, u) call add_array_r2(cv3a_compress_list, v1, v) call add_array_r2(cv3a_compress_list, gz1, gz) call add_array_r2(cv3a_compress_list, h1, h) call add_array_r2(cv3a_compress_list, th1_wake, th_wake) call add_array_r2(cv3a_compress_list, lv1, lv) call add_array_r2(cv3a_compress_list, lf1, lf) call add_array_r2(cv3a_compress_list, cpn1, cpn) call add_array_r2(cv3a_compress_list, p1, p) call add_array_r2(cv3a_compress_list, ph1, ph) call add_array_r2(cv3a_compress_list, tv1, tv) call add_array_r2(cv3a_compress_list, tp1, tp) call add_array_r2(cv3a_compress_list, tvp1, tvp) call add_array_r2(cv3a_compress_list, clw1, clw) call add_array_r2(cv3a_compress_list, h1_wake, h_wake) call add_array_r2(cv3a_compress_list, lv1_wake, lv_wake) call add_array_r2(cv3a_compress_list, lf1_wake, lf_wake) call add_array_r2(cv3a_compress_list, cpn1_wake, cpn_wake) call add_array_r2(cv3a_compress_list, tv1_wake, tv_wake) call add_array_r2(cv3a_compress_list, sig1, sig) call add_array_r1(cv3a_compress_list, sig1(:, nd), sig(:, nd)) call add_array_r2(cv3a_compress_list, w01, w0) call add_array_r1(cv3a_compress_list, Ale1, Ale) call add_array_r1(cv3a_compress_list, Alp1, Alp) call add_array_r2(cv3a_compress_list, omega1, omega) call cv3a_compress(len, (iflag1 == 0), cv3a_compress_list, compress_data) ncum = compress_data%ncum call exit_profile("cv3a_compress") call cv3a_driver_compressed(nloc, nd, ntra, & debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, & ! Input fields delt, tau_cld_cv, coefw_cld_cv, nk, icbs, tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp, wghti, th, t, q, qs, t_wake, q_wake, qs_wake, u, v, gz, h, th_wake, lv, lf, cpn, p, tv, tp, h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, omega, ph, & ! Input/output iflag, icb, plcl, sig, w0, tvp, clw, & ! Output fields inb, precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag, ft, fq, fu, fv, ma, mip, upwd, dnwd, dnwd0, qcondc, ftd, fqd, asupmax, da, mp, d1a, dam, qta, evap, ep, eplaMm, wdtrainA, wdtrainS, wdtrainM, qtc, sigt, vprecip, vprecipi, phi, phi2, sigij, elij, epmlmMm) ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ ! --- UNCOMPRESS THE FIELDS ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ call enter_profile("cv3a_uncompress") call driver_log('cv3a_uncompress') call add_array_i1(cv3a_uncompress_list, iflag, iflag1, init=.false.) call add_array_i1(cv3a_uncompress_list, icb, kbas1) call add_array_i1(cv3a_uncompress_list, inb, ktop1) call add_array_r1(cv3a_uncompress_list, precip, precip1) call add_array_r1(cv3a_uncompress_list, cbmf, cbmf1) call add_array_r1(cv3a_uncompress_list, plcl, plcl1, init=.false.) call add_array_r1(cv3a_uncompress_list, plfc, plfc1) call add_array_r1(cv3a_uncompress_list, wbeff, wbeff1) call add_array_r2(cv3a_uncompress_list, sig, sig1) call add_array_r2(cv3a_uncompress_list, w0, w01) call add_array_r1(cv3a_uncompress_list, ptop2, ptop21) call add_array_r2(cv3a_uncompress_list, ft, ft1) call add_array_r2(cv3a_uncompress_list, fq, fq1) call add_array_r2(cv3a_uncompress_list, fu, fu1) call add_array_r2(cv3a_uncompress_list, fv, fv1) call add_array_r1(cv3a_uncompress_list, sigd, sigd1) call add_array_r2(cv3a_uncompress_list, ma, ma1) call add_array_r2(cv3a_uncompress_list, mip, mip1) call add_array_r2(cv3a_uncompress_list, vprecip, vprecip1) call add_array_r2(cv3a_uncompress_list, vprecipi, vprecipi1) call add_array_r2(cv3a_uncompress_list, upwd, upwd1) call add_array_r2(cv3a_uncompress_list, dnwd, dnwd1) call add_array_r2(cv3a_uncompress_list, dnwd0, dnwd01) call add_array_r2(cv3a_uncompress_list, qcondc, qcondc1) call add_array_r1(cv3a_uncompress_list, wd, wd1) cape1(:) = Cape_noconv call add_array_r1(cv3a_uncompress_list, cape, cape1, init=.false.) cin1(:) = Cin_noconv call add_array_r1(cv3a_uncompress_list, cin, cin1, init=.false.) call add_array_r2(cv3a_uncompress_list, tvp, tvp1, init=.false.) call add_array_r2(cv3a_uncompress_list, ftd, ftd1) call add_array_r2(cv3a_uncompress_list, fqd, fqd1) call add_array_r1(cv3a_uncompress_list, Plim1, Plim11) call add_array_r1(cv3a_uncompress_list, plim2, plim21) call add_array_r2(cv3a_uncompress_list, asupmax, asupmax1) call add_array_r1(cv3a_uncompress_list, supmax0, supmax01) call add_array_r1(cv3a_uncompress_list, asupmaxmin, asupmaxmin1) call add_array_r2(cv3a_uncompress_list, da, da1) call add_array_r3(cv3a_uncompress_list, phi, phi1) call add_array_r2(cv3a_uncompress_list, mp, mp1) call add_array_r3(cv3a_uncompress_list, phi2, phi21) call add_array_r2(cv3a_uncompress_list, d1a, d1a1) call add_array_r2(cv3a_uncompress_list, dam, dam1) call add_array_r3(cv3a_uncompress_list, sigij, sigij1) call add_array_r2(cv3a_uncompress_list, qta, qta1) call add_array_r2(cv3a_uncompress_list, clw, clw1, init=.false.) call add_array_r3(cv3a_uncompress_list, elij, elij1) call add_array_r2(cv3a_uncompress_list, evap, evap1) call add_array_r2(cv3a_uncompress_list, ep, ep1) call add_array_r3(cv3a_uncompress_list, epmlmMm, epmlmMm1) call add_array_r2(cv3a_uncompress_list, eplaMm, eplaMm1) call add_array_r2(cv3a_uncompress_list, wdtrainA, wdtrainA1) call add_array_r2(cv3a_uncompress_list, wdtrainS, wdtrainS1) call add_array_r2(cv3a_uncompress_list, wdtrainM, wdtrainM1) call add_array_r2(cv3a_uncompress_list, qtc, qtc1) call add_array_r2(cv3a_uncompress_list, sigt, sigt1) call add_array_r1(cv3a_uncompress_list, epmax_diag, epmax_diag1) call add_array_r1(cv3a_uncompress_list, sig(:, nd), sig1(:, nd)) call cv3a_uncompress(len, compress_data, cv3a_uncompress_list) call exit_profile("cv3a_uncompress") end block else !compress nloc = len coef_convective(:) = merge(1, 0, iflag1(:) == 0) call cv3a_driver_compressed(nloc, nd, ntra, & debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, & ! Input fields delt, tau_cld_cv, coefw_cld_cv, nk1, icbs1, tnk1, qnk1, gznk1, hnk1, unk1, vnk1, pbase1, buoybase1, s1_wake, Ale1, Alp1, wghti1, th1, t1, q1, qs1, t1_wake, q1_wake, qs1_wake, u1, v1, gz1, h1, th1_wake, lv1, lf1, cpn1, p1, tv1, tp1, h1_wake, lv1_wake, lf1_wake, cpn1_wake, tv1_wake, omega1, ph1, & ! Input/output iflag1, kbas1, plcl1, sig1, w01, tvp1, clw1, & ! Output fields ktop1, precip1, cbmf1, plfc1, wbeff1, ptop21, sigd1, wd1, cape1, cin1, Plim11, plim21, supmax01, asupmaxmin1, epmax_diag1, ft1, fq1, fu1, fv1, ma1, mip1, upwd1, dnwd1, dnwd01, qcondc1, ftd1, fqd1, asupmax1, da1, mp1, d1a1, dam1, qta1, evap1, ep1, eplaMm1, wdtrainA1, wdtrainS1, wdtrainM1, qtc1, sigt1, vprecip1, vprecipi1, phi1, phi21, sigij1, elij1, epmlmMm1) ! Reset Cin1 and Cape1 to their default values for non-convective grid-points. Cin1(:) = Cin1(:)*coef_convective(:) + Cin_noconv*(1.-coef_convective(:)) Cape1(:) = Cape1(:)*coef_convective(:) + Cape_noconv*(1.-coef_convective(:)) precip1(:) = precip1(:)*coef_convective(:) kbas1(:) = kbas1(:)*coef_convective(:) ktop1(:) = ktop1(:)*coef_convective(:) sigd1(:) = sigd1(:)*coef_convective(:) wbeff1(:) = wbeff1(:)*coef_convective(:) DO k = 1, nd qcondc1(:, k) = qcondc1(:, k)*coef_convective(:) vprecip1(:, k) = vprecip1(:, k)*coef_convective(:) vprecipi1(:, k) = vprecipi1(:, k)*coef_convective(:) ENDDO endif IF (prt_level >= 10) THEN Print *, 'cva_driver after cv3_uncompress:ft1(1) , ftd1(1) ', & ft1(igout, 1), ftd1(igout, 1) Print *, 'cva_driver after cv3_uncompress:fq1(1) , fqd1(1) ', & fq1(igout, 1), fqd1(igout, 1) ENDIF IF (debut) THEN PRINT *, ' cv_uncompress -> ' debut = .FALSE. END IF ftra1(:, :, :) = 0 call exit_profile("cv3a_driver") END SUBROUTINE subroutine cv3a_driver_compressed(nloc, nd, ntra, & debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, & ! Input fields delt, tau_cld_cv, coefw_cld_cv, nk, icbs, tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp, wghti, th, t, q, qs, t_wake, q_wake, qs_wake, u, v, gz, h, th_wake, lv, lf, cpn, p, tv, tp, h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, omega, ph, & ! Input/output iflag, icb, plcl, sig, w0, tvp, clw, & ! Output fields inb, precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag, ft, fq, fu, fv, ma, mip, upwd, dnwd, dnwd0, qcondc, ftd, fqd, asupmax, da, mp, d1a, dam, qta, evap, ep, eplaMm, wdtrainA, wdtrainS, wdtrainM, qtc, sigt, vprecip, vprecipi, phi, phi2, sigij, elij, epmlmMm) use profiling_physic_mod USE cv3p_mixing_mod, ONLY : cv3p_mixing integer, intent(in) :: nloc, nd, ntra ! Arrays sizes logical, intent(in) :: debut, ok_conserv_q integer, intent(in) :: iflag_mix, iflag_clos, prt_level real, intent(in) :: delt, tau_cld_cv, coefw_cld_cv ! Input arrays integer, dimension(nloc), intent(in) :: nk, icbs real, dimension(nloc), intent(in) :: tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp real, dimension(nloc, nd), intent(in) :: wghti, th, t, q, qs, t_wake, q_wake, qs_wake, u, v, gz, h, th_wake, lv, lf, cpn, p, tv, tp, h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, omega real, dimension(nloc, nd + 1), intent(in) :: ph ! Input/output arrays integer, dimension(nloc), intent(inout) :: iflag, icb real, dimension(nloc), intent(inout) :: plcl real, dimension(nloc, nd), intent(inout) :: sig, w0, tvp, clw ! Output arrays integer, dimension(nloc), intent(out) :: inb real, dimension(nloc), intent(out) :: precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag real, dimension(nloc, nd), intent(out) :: ft, fq, fu, fv, ma, mip, upwd, dnwd, dnwd0, qcondc, ftd, fqd, asupmax, da, mp, d1a, dam, qta, evap, ep, eplaMm, wdtrainA, wdtrainS, wdtrainM, qtc, sigt real, dimension(nloc, nd + 1), intent(out) :: vprecip, vprecipi real, dimension(nloc, nd, nd), intent(out) :: phi, phi2, sigij, elij, epmlmMm ! Local fields INTEGER nent(nloc, nd) REAL hp(nloc, nd), sigp(nloc, nd), qpreca(nloc, nd) REAL buoy(nloc, nd) REAL m(nloc, nd) REAL mm(nloc, nd) REAL ment(nloc, nd, nd) REAL qent(nloc, nd, nd) REAL hent(nloc, nd, nd) REAL uent(nloc, nd, nd), vent(nloc, nd, nd) REAL ments(nloc, nd, nd), qents(nloc, nd, nd) REAL supmax(nloc, nd) REAL coef_clos(nloc) REAL, DIMENSION(nloc, nd) :: qp, up, vp REAL, DIMENSION(nloc, nd) :: wt, water REAL, DIMENSION(nloc, nd) :: ice, fondue, b REAL, DIMENSION(nloc, nd) :: frac_a, frac_s, faci REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra) REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra) LOGICAL ok_inhib ! True => possible inhibition of convection by dryness integer :: k integer :: ncum ! Number of convective cells logical, save :: timers_first = .true. !$omp threadprivate(timers_first) INTEGER, PARAMETER :: igout = 1 ncum = nloc call enter_profile("cv3a_compresd") IF (ncum > 0) THEN ! ------------------------------------------------------------------- ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part : ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES ! --- & ! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE ! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD ! --- & ! --- FIND THE LEVEL OF NEUTRAL BUOYANCY ! ------------------------------------------------------------------- call driver_log('cv3_undilute2') qta(:,:) = 0 ep(:,:) = 0 CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, & tnk, qnk, gznk, hnk, t, q, qs, gz, & p, ph, h, tv, lv, lf, pbase, buoybase, plcl, & inb, tp, tvp, clw, hp, ep, sigp, buoy, & frac_a, frac_s, qpreca, qta) ! epmax_cape ! on recalcule ep et hp call driver_log('cv3_epmax_cape') epmax_diag(:) = 0 call cv3_epmax_fn_cape(nloc, ncum, nd & , ep, hp, icb, inb, clw, nk, t, h, hnk, lv, lf, frac_s & , pbase, p, ph, tv, buoy, sig, w0, iflag & , epmax_diag) ! ------------------------------------------------------------------- ! --- MIXING(1) (if iflag_mix .ge. 1) ! ------------------------------------------------------------------- elij(:,:,:) = 0 supmax(:,:) = 0 IF (iflag_mix >= 1) THEN call driver_log('cv3p_mixing') CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, & ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qta, & unk, vnk, hp, tv, tvp, ep, clw, sig, & ment, qent, hent, uent, vent, nent, & sigij, elij, supmax, ments, qents, traent) END IF ! ------------------------------------------------------------------- ! --- CLOSURE ! ------------------------------------------------------------------- cape(:) = -1 cin(:) = -100000. ptop2(:) = 0 coef_clos(:) = 1. Plim1(:) = 0 plim2(:) = 0 supmax0(:) = 0 asupmaxmin(:) = 0 cbmf(:) = 0 plfc(:) = 0 wbeff(:) = 0 asupmax(:,:) = 0 ok_inhib = (iflag_mix == 2) IF (iflag_clos == 0) THEN call driver_log('cv3_closure') CALL cv3_closure(nloc, ncum, nd, icb, inb, & pbase, p, ph, tv, buoy, & sig, w0, cape, m, iflag) END IF ! iflag_clos==0 IF (iflag_clos == 1) PRINT *, ' pas d appel cv3p_closure' IF (iflag_clos == 2) THEN call driver_log('cv3p1_closure') CALL cv3p1_closure(nloc, ncum, nd, icb, inb, & pbase, plcl, p, ph, tv, tvp, buoy, & supmax, ok_inhib, Ale, Alp, omega, & sig, w0, ptop2, cape, cin, m, iflag, coef_clos, & Plim1, plim2, asupmax, supmax0, & asupmaxmin, cbmf, plfc, wbeff) if (prt_level >= 10) PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1) END IF ! iflag_clos==2 IF (iflag_clos == 3) THEN call driver_log('cv3p2_closure') CALL cv3p2_closure(nloc, ncum, nd, icb, inb, & pbase, plcl, p, ph, tv, tvp, buoy, & supmax, ok_inhib, Ale, Alp, omega, & sig, w0, ptop2, cape, cin, m, iflag, coef_clos, & Plim1, plim2, asupmax, supmax0, & asupmaxmin, cbmf, plfc, wbeff) if (prt_level >= 10) PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1) END IF ! iflag_clos==3 ! ------------------------------------------------------------------- ! --- MIXING(2) ! ------------------------------------------------------------------- IF (iflag_mix == 0) THEN call driver_log('cv3_mixing') CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, & ! na->nd ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qnk, & unk, vnk, hp, tv, tvp, ep, clw, m, sig, & ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent) CALL zilch(hent, nloc*nd*nd) ELSE mm(:, :) = m(:, :) CALL cv3_mixscale(nloc, ncum, nd, ment, mm) IF (debut) PRINT *, ' cv3_mixscale-> ' END IF IF (debut) PRINT *, ' cv_mixing ->' ! ------------------------------------------------------------------- ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS ! ------------------------------------------------------------------- IF (debut) PRINT *, ' cva_driver -> cv3_unsat ' call driver_log('cv3_unsat') CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, & t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, & th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, & ep, sigp, clw, frac_s, qpreca, frac_a, qta, & m, ment, elij, delt, plcl, coef_clos, & mp, qp, up, vp, trap, wt, water, evap, fondue, ice, & faci, b, sigd, & wdtrainA, wdtrainS, wdtrainM) IF (prt_level >= 10) THEN Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue ' DO k = 1, nd write (6, '(i4,5(1x,e13.6))'), & k, mp(igout, k), water(igout, k), ice(igout, k), & evap(igout, k), fondue(igout, k) ENDDO Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainS, wdtrainM ' DO k = 1, nd write (6, '(i4,3(1x,e13.6))'), & k, wdtrainA(igout, k), wdtrainS(igout, k), wdtrainM(igout, k) ENDDO ENDIF IF (debut) PRINT *, 'cv_unsat-> ' ! ------------------------------------------------------------------- ! YIELD ! (tendencies, precipitation, variables of interface with other processes, etc) ! ------------------------------------------------------------------- call driver_log('cv3_yield') CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, & icb, inb, delt, & t, q, t_wake, q_wake, s_wake, u, v, tra, & gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, & ep, clw, qpreca, m, tp, mp, qp, up, vp, trap, & wt, water, ice, evap, fondue, faci, b, sigd, & ment, qent, hent, iflag_mix, uent, vent, & nent, elij, traent, sig, & tv, tvp, wghti, & iflag, precip, Vprecip, Vprecipi, ft, fq, fu, fv, ftra, & cbmf, upwd, dnwd, dnwd0, ma, mip, & qcondc, wd, & ftd, fqd, qta, qtc, sigt, tau_cld_cv, coefw_cld_cv) ! Test conseravtion de l'eau IF (debut) PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1) IF (prt_level >= 10) THEN Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', & ft(igout, 1), ftd(igout, 1) Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', & fq(igout, 1), fqd(igout, 1) ENDIF !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ !--- passive tracers !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ call driver_log('cv3_tracer') sigij(:,:,:) = 0 vprecip(:,:) = 0 ! GLITCHY : vprecip is unused in cv3_tracer and sigij is intent in CALL cv3_tracer(nloc, -1, ncum, nd, nd, & ment, sigij, da, phi, phi2, d1a, dam, & ep, vprecip, elij, clw, epmlmMm, eplaMm, & icb, inb) vprecipi = 0 END IF ! ncum>0 call exit_profile("cv3a_compresd") end subroutine subroutine driver_log(message) use print_control_mod, only: prt_level character(*) :: message if (prt_level >= 9) PRINT *, 'cva_driver ->', message end subroutine end module