[524] | 1 | ! |
---|
| 2 | ! $Header$ |
---|
| 3 | ! |
---|
| 4 | SUBROUTINE cv_driver(len,nd,ndp1,ntra,iflag_con, |
---|
| 5 | & t1,q1,qs1,u1,v1,tra1, |
---|
| 6 | & p1,ph1,iflag1,ft1,fq1,fu1,fv1,ftra1, |
---|
[619] | 7 | & precip1,VPrecip1, |
---|
[524] | 8 | & cbmf1,sig1,w01, |
---|
[619] | 9 | & icb1,inb1, |
---|
| 10 | & delt,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1, |
---|
[1774] | 11 | & da1,phi1,mp1,phi21,d1a1,dam1,sij1,clw1,elij1, ! RomP |
---|
| 12 | & evap1,ep1,epmlmMm1,eplaMm1, ! RomP |
---|
| 13 | & wdtrainA1,wdtrainM1) ! RomP |
---|
[524] | 14 | C |
---|
[766] | 15 | USE dimphy |
---|
[524] | 16 | implicit none |
---|
| 17 | C |
---|
| 18 | C.............................START PROLOGUE............................ |
---|
| 19 | C |
---|
[1774] | 20 | ! |
---|
| 21 | ! All argument names (except len,nd,ntra,nloc,delt and the flags) have a "1" appended. |
---|
| 22 | ! The "1" is removed for the corresponding compressed (local) variables. |
---|
| 23 | ! |
---|
[524] | 24 | C PARAMETERS: |
---|
| 25 | C Name Type Usage Description |
---|
| 26 | C ---------- ---------- ------- ---------------------------- |
---|
| 27 | C |
---|
| 28 | C len Integer Input first (i) dimension |
---|
| 29 | C nd Integer Input vertical (k) dimension |
---|
| 30 | C ndp1 Integer Input nd + 1 |
---|
| 31 | C ntra Integer Input number of tracors |
---|
| 32 | C iflag_con Integer Input version of convect (3/4) |
---|
| 33 | C t1 Real Input temperature |
---|
| 34 | C q1 Real Input specific hum |
---|
| 35 | C qs1 Real Input sat specific hum |
---|
| 36 | C u1 Real Input u-wind |
---|
| 37 | C v1 Real Input v-wind |
---|
| 38 | C tra1 Real Input tracors |
---|
| 39 | C p1 Real Input full level pressure |
---|
| 40 | C ph1 Real Input half level pressure |
---|
| 41 | C iflag1 Integer Output flag for Emanuel conditions |
---|
| 42 | C ft1 Real Output temp tend |
---|
| 43 | C fq1 Real Output spec hum tend |
---|
| 44 | C fu1 Real Output u-wind tend |
---|
| 45 | C fv1 Real Output v-wind tend |
---|
| 46 | C ftra1 Real Output tracor tend |
---|
| 47 | C precip1 Real Output precipitation |
---|
[619] | 48 | C VPrecip1 Real Output vertical profile of precipitations |
---|
[524] | 49 | C cbmf1 Real Output cloud base mass flux |
---|
| 50 | C sig1 Real In/Out section adiabatic updraft |
---|
| 51 | C w01 Real In/Out vertical velocity within adiab updraft |
---|
| 52 | C delt Real Input time step |
---|
| 53 | C Ma1 Real Output mass flux adiabatic updraft |
---|
| 54 | C upwd1 Real Output total upward mass flux (adiab+mixed) |
---|
| 55 | C dnwd1 Real Output saturated downward mass flux (mixed) |
---|
| 56 | C dnwd01 Real Output unsaturated downward mass flux |
---|
| 57 | C qcondc1 Real Output in-cld mixing ratio of condensed water |
---|
| 58 | C wd1 Real Output downdraft velocity scale for sfc fluxes |
---|
| 59 | C cape1 Real Output CAPE |
---|
[1774] | 60 | ! |
---|
| 61 | ! wdtrainA1 Real Output precipitation detrained from adiabatic draught; |
---|
| 62 | ! used in tracer transport (cvltr) |
---|
| 63 | ! wdtrainM1 Real Output precipitation detrained from mixed draughts; |
---|
| 64 | ! used in tracer transport (cvltr) |
---|
| 65 | ! da1 Real Output used in tracer transport (cvltr) |
---|
| 66 | ! phi1 Real Output used in tracer transport (cvltr) |
---|
| 67 | ! mp1 Real Output used in tracer transport (cvltr) |
---|
| 68 | ! |
---|
| 69 | ! phi21 Real Output used in tracer transport (cvltr) |
---|
| 70 | ! |
---|
| 71 | ! d1a1 Real Output used in tracer transport (cvltr) |
---|
| 72 | ! dam1 Real Output used in tracer transport (cvltr) |
---|
| 73 | ! |
---|
| 74 | ! evap1 Real Output |
---|
| 75 | ! ep1 Real Output |
---|
| 76 | ! sij1 Real Output |
---|
| 77 | ! elij1 Real Output |
---|
[524] | 78 | C |
---|
| 79 | C S. Bony, Mar 2002: |
---|
| 80 | C * Several modules corresponding to different physical processes |
---|
| 81 | C * Several versions of convect may be used: |
---|
| 82 | C - iflag_con=3: version lmd (previously named convect3) |
---|
| 83 | C - iflag_con=4: version 4.3b (vect. version, previously convect1/2) |
---|
| 84 | C + tard: - iflag_con=5: version lmd with ice (previously named convectg) |
---|
| 85 | C S. Bony, Oct 2002: |
---|
| 86 | C * Vectorization of convect3 (ie version lmd) |
---|
| 87 | C |
---|
| 88 | C..............................END PROLOGUE............................. |
---|
| 89 | c |
---|
| 90 | c |
---|
[766] | 91 | cym#include "dimensions.h" |
---|
| 92 | cym#include "dimphy.h" |
---|
[1742] | 93 | c |
---|
| 94 | c Input |
---|
[524] | 95 | integer len |
---|
| 96 | integer nd |
---|
| 97 | integer ndp1 |
---|
| 98 | integer noff |
---|
| 99 | integer iflag_con |
---|
| 100 | integer ntra |
---|
[1742] | 101 | real delt |
---|
[524] | 102 | real t1(len,nd) |
---|
| 103 | real q1(len,nd) |
---|
| 104 | real qs1(len,nd) |
---|
| 105 | real u1(len,nd) |
---|
| 106 | real v1(len,nd) |
---|
[1742] | 107 | real tra1(len,nd,ntra) |
---|
[524] | 108 | real p1(len,nd) |
---|
| 109 | real ph1(len,ndp1) |
---|
[1742] | 110 | c |
---|
| 111 | c Output |
---|
[524] | 112 | integer iflag1(len) |
---|
| 113 | real ft1(len,nd) |
---|
| 114 | real fq1(len,nd) |
---|
| 115 | real fu1(len,nd) |
---|
| 116 | real fv1(len,nd) |
---|
[1742] | 117 | real ftra1(len,nd,ntra) |
---|
[524] | 118 | real precip1(len) |
---|
| 119 | real cbmf1(len) |
---|
[1742] | 120 | real sig1(klon,klev) |
---|
| 121 | real w01(klon,klev) |
---|
| 122 | real VPrecip1(len,nd+1) |
---|
| 123 | real evap1(len,nd) !RomP |
---|
| 124 | real ep1(len,nd) !RomP |
---|
[524] | 125 | real Ma1(len,nd) |
---|
| 126 | real upwd1(len,nd) |
---|
| 127 | real dnwd1(len,nd) |
---|
| 128 | real dnwd01(len,nd) |
---|
| 129 | |
---|
| 130 | real qcondc1(len,nd) ! cld |
---|
| 131 | real wd1(len) ! gust |
---|
| 132 | real cape1(len) |
---|
| 133 | |
---|
[1742] | 134 | ! RomP >>> |
---|
| 135 | real wdtrainA1(len,nd), wdtrainM1(len,nd) |
---|
| 136 | real sij1(len,nd,nd),elij1(len,nd,nd) |
---|
[619] | 137 | real da1(len,nd),phi1(len,nd,nd),mp1(len,nd) |
---|
[524] | 138 | |
---|
[1742] | 139 | real phi21(len,nd,nd) |
---|
| 140 | real d1a1(len,nd), dam1(len,nd) |
---|
| 141 | real epmlmMm1(len,nd,nd),eplaMm1(len,nd) |
---|
| 142 | ! RomP <<< |
---|
[1774] | 143 | ! |
---|
[524] | 144 | !------------------------------------------------------------------- |
---|
[1774] | 145 | ! Original Prologue by Kerry Emanuel. |
---|
| 146 | !------------------------------------------------------------------- |
---|
[524] | 147 | ! --- ARGUMENTS |
---|
| 148 | !------------------------------------------------------------------- |
---|
| 149 | ! --- On input: |
---|
| 150 | ! |
---|
| 151 | ! t: Array of absolute temperature (K) of dimension ND, with first |
---|
| 152 | ! index corresponding to lowest model level. Note that this array |
---|
| 153 | ! will be altered by the subroutine if dry convective adjustment |
---|
| 154 | ! occurs and if IPBL is not equal to 0. |
---|
| 155 | ! |
---|
| 156 | ! q: Array of specific humidity (gm/gm) of dimension ND, with first |
---|
| 157 | ! index corresponding to lowest model level. Must be defined |
---|
| 158 | ! at same grid levels as T. Note that this array will be altered |
---|
| 159 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
| 160 | ! |
---|
| 161 | ! qs: Array of saturation specific humidity of dimension ND, with first |
---|
| 162 | ! index corresponding to lowest model level. Must be defined |
---|
| 163 | ! at same grid levels as T. Note that this array will be altered |
---|
| 164 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
| 165 | ! |
---|
| 166 | ! u: Array of zonal wind velocity (m/s) of dimension ND, witth first |
---|
| 167 | ! index corresponding with the lowest model level. Defined at |
---|
| 168 | ! same levels as T. Note that this array will be altered if |
---|
| 169 | ! dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
| 170 | ! |
---|
| 171 | ! v: Same as u but for meridional velocity. |
---|
| 172 | ! |
---|
| 173 | ! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA), |
---|
| 174 | ! where NTRA is the number of different tracers. If no |
---|
| 175 | ! convective tracer transport is needed, define a dummy |
---|
| 176 | ! input array of dimension (ND,1). Tracers are defined at |
---|
| 177 | ! same vertical levels as T. Note that this array will be altered |
---|
| 178 | ! if dry convective adjustment occurs and if IPBL is not equal to 0. |
---|
| 179 | ! |
---|
| 180 | ! p: Array of pressure (mb) of dimension ND, with first |
---|
| 181 | ! index corresponding to lowest model level. Must be defined |
---|
| 182 | ! at same grid levels as T. |
---|
| 183 | ! |
---|
| 184 | ! ph: Array of pressure (mb) of dimension ND+1, with first index |
---|
| 185 | ! corresponding to lowest level. These pressures are defined at |
---|
| 186 | ! levels intermediate between those of P, T, Q and QS. The first |
---|
| 187 | ! value of PH should be greater than (i.e. at a lower level than) |
---|
| 188 | ! the first value of the array P. |
---|
| 189 | ! |
---|
| 190 | ! nl: The maximum number of levels to which convection can penetrate, plus 1. |
---|
| 191 | ! NL MUST be less than or equal to ND-1. |
---|
| 192 | ! |
---|
| 193 | ! delt: The model time step (sec) between calls to CONVECT |
---|
| 194 | ! |
---|
| 195 | !---------------------------------------------------------------------------- |
---|
| 196 | ! --- On Output: |
---|
| 197 | ! |
---|
| 198 | ! iflag: An output integer whose value denotes the following: |
---|
| 199 | ! VALUE INTERPRETATION |
---|
| 200 | ! ----- -------------- |
---|
| 201 | ! 0 Moist convection occurs. |
---|
| 202 | ! 1 Moist convection occurs, but a CFL condition |
---|
| 203 | ! on the subsidence warming is violated. This |
---|
| 204 | ! does not cause the scheme to terminate. |
---|
| 205 | ! 2 Moist convection, but no precip because ep(inb) lt 0.0001 |
---|
| 206 | ! 3 No moist convection because new cbmf is 0 and old cbmf is 0. |
---|
| 207 | ! 4 No moist convection; atmosphere is not |
---|
| 208 | ! unstable |
---|
| 209 | ! 6 No moist convection because ihmin le minorig. |
---|
| 210 | ! 7 No moist convection because unreasonable |
---|
| 211 | ! parcel level temperature or specific humidity. |
---|
| 212 | ! 8 No moist convection: lifted condensation |
---|
| 213 | ! level is above the 200 mb level. |
---|
| 214 | ! 9 No moist convection: cloud base is higher |
---|
| 215 | ! then the level NL-1. |
---|
| 216 | ! |
---|
| 217 | ! ft: Array of temperature tendency (K/s) of dimension ND, defined at same |
---|
| 218 | ! grid levels as T, Q, QS and P. |
---|
| 219 | ! |
---|
| 220 | ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND, |
---|
| 221 | ! defined at same grid levels as T, Q, QS and P. |
---|
| 222 | ! |
---|
| 223 | ! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND, |
---|
| 224 | ! defined at same grid levels as T. |
---|
| 225 | ! |
---|
| 226 | ! fv: Same as FU, but for forcing of meridional velocity. |
---|
| 227 | ! |
---|
| 228 | ! ftra: Array of forcing of tracer content, in tracer mixing ratio per |
---|
| 229 | ! second, defined at same levels as T. Dimensioned (ND,NTRA). |
---|
| 230 | ! |
---|
| 231 | ! precip: Scalar convective precipitation rate (mm/day). |
---|
| 232 | ! |
---|
[619] | 233 | ! VPrecip: Vertical profile of convective precipitation (kg/m2/s). |
---|
| 234 | ! |
---|
[524] | 235 | ! wd: A convective downdraft velocity scale. For use in surface |
---|
| 236 | ! flux parameterizations. See convect.ps file for details. |
---|
| 237 | ! |
---|
| 238 | ! tprime: A convective downdraft temperature perturbation scale (K). |
---|
| 239 | ! For use in surface flux parameterizations. See convect.ps |
---|
| 240 | ! file for details. |
---|
| 241 | ! |
---|
| 242 | ! qprime: A convective downdraft specific humidity |
---|
| 243 | ! perturbation scale (gm/gm). |
---|
| 244 | ! For use in surface flux parameterizations. See convect.ps |
---|
| 245 | ! file for details. |
---|
| 246 | ! |
---|
| 247 | ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST |
---|
| 248 | ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT |
---|
| 249 | ! ITS NEXT CALL. That is, the value of CBMF must be "remembered" |
---|
| 250 | ! by the calling program between calls to CONVECT. |
---|
| 251 | ! |
---|
| 252 | ! det: Array of detrainment mass flux of dimension ND. |
---|
| 253 | ! |
---|
| 254 | !------------------------------------------------------------------- |
---|
| 255 | c |
---|
| 256 | c Local arrays |
---|
| 257 | c |
---|
| 258 | |
---|
| 259 | integer i,k,n,il,j |
---|
| 260 | integer icbmax |
---|
| 261 | integer nk1(klon) |
---|
| 262 | integer icb1(klon) |
---|
[619] | 263 | integer inb1(klon) |
---|
[524] | 264 | integer icbs1(klon) |
---|
| 265 | |
---|
| 266 | real plcl1(klon) |
---|
| 267 | real tnk1(klon) |
---|
| 268 | real qnk1(klon) |
---|
| 269 | real gznk1(klon) |
---|
| 270 | real pnk1(klon) |
---|
| 271 | real qsnk1(klon) |
---|
| 272 | real pbase1(klon) |
---|
| 273 | real buoybase1(klon) |
---|
| 274 | |
---|
| 275 | real lv1(klon,klev) |
---|
| 276 | real cpn1(klon,klev) |
---|
| 277 | real tv1(klon,klev) |
---|
| 278 | real gz1(klon,klev) |
---|
| 279 | real hm1(klon,klev) |
---|
| 280 | real h1(klon,klev) |
---|
| 281 | real tp1(klon,klev) |
---|
| 282 | real tvp1(klon,klev) |
---|
| 283 | real clw1(klon,klev) |
---|
| 284 | real th1(klon,klev) |
---|
| 285 | c |
---|
| 286 | integer ncum |
---|
| 287 | c |
---|
| 288 | c (local) compressed fields: |
---|
| 289 | c |
---|
[766] | 290 | cym integer nloc |
---|
| 291 | cym parameter (nloc=klon) ! pour l'instant |
---|
| 292 | #define nloc klon |
---|
[524] | 293 | integer idcum(nloc) |
---|
| 294 | integer iflag(nloc),nk(nloc),icb(nloc) |
---|
| 295 | integer nent(nloc,klev) |
---|
| 296 | integer icbs(nloc) |
---|
| 297 | integer inb(nloc), inbis(nloc) |
---|
| 298 | |
---|
| 299 | real cbmf(nloc),plcl(nloc),tnk(nloc),qnk(nloc),gznk(nloc) |
---|
| 300 | real t(nloc,klev),q(nloc,klev),qs(nloc,klev) |
---|
| 301 | real u(nloc,klev),v(nloc,klev) |
---|
| 302 | real gz(nloc,klev),h(nloc,klev),lv(nloc,klev),cpn(nloc,klev) |
---|
| 303 | real p(nloc,klev),ph(nloc,klev+1),tv(nloc,klev),tp(nloc,klev) |
---|
| 304 | real clw(nloc,klev) |
---|
| 305 | real dph(nloc,klev) |
---|
| 306 | real pbase(nloc), buoybase(nloc), th(nloc,klev) |
---|
| 307 | real tvp(nloc,klev) |
---|
| 308 | real sig(nloc,klev), w0(nloc,klev) |
---|
| 309 | real hp(nloc,klev), ep(nloc,klev), sigp(nloc,klev) |
---|
| 310 | real frac(nloc), buoy(nloc,klev) |
---|
| 311 | real cape(nloc) |
---|
| 312 | real m(nloc,klev), ment(nloc,klev,klev), qent(nloc,klev,klev) |
---|
| 313 | real uent(nloc,klev,klev), vent(nloc,klev,klev) |
---|
| 314 | real ments(nloc,klev,klev), qents(nloc,klev,klev) |
---|
| 315 | real sij(nloc,klev,klev), elij(nloc,klev,klev) |
---|
[619] | 316 | real qp(nloc,klev), up(nloc,klev), vp(nloc,klev) |
---|
[524] | 317 | real wt(nloc,klev), water(nloc,klev), evap(nloc,klev) |
---|
| 318 | real b(nloc,klev), ft(nloc,klev), fq(nloc,klev) |
---|
| 319 | real fu(nloc,klev), fv(nloc,klev) |
---|
| 320 | real upwd(nloc,klev), dnwd(nloc,klev), dnwd0(nloc,klev) |
---|
| 321 | real Ma(nloc,klev), mike(nloc,klev), tls(nloc,klev) |
---|
| 322 | real tps(nloc,klev), qprime(nloc), tprime(nloc) |
---|
| 323 | real precip(nloc) |
---|
[619] | 324 | real VPrecip(nloc,klev+1) |
---|
[524] | 325 | real tra(nloc,klev,ntra), trap(nloc,klev,ntra) |
---|
| 326 | real ftra(nloc,klev,ntra), traent(nloc,klev,klev,ntra) |
---|
| 327 | real qcondc(nloc,klev) ! cld |
---|
| 328 | real wd(nloc) ! gust |
---|
[1774] | 329 | ! |
---|
| 330 | ! RomP >>> |
---|
| 331 | real da(nloc,klev),phi(nloc,klev,klev),mp(nloc,klev) |
---|
| 332 | real epmlmMm(nloc,klev,klev),eplaMm(nloc,klev) |
---|
| 333 | real phi2(nloc,klev,klev) |
---|
| 334 | real d1a(nloc,klev), dam(nloc,klev) |
---|
| 335 | real wdtrainA(nloc,klev),wdtrainM(nloc,klev) |
---|
| 336 | real sigd(nloc) |
---|
| 337 | ! RomP <<< |
---|
[524] | 338 | |
---|
[766] | 339 | nent(:,:)=0 |
---|
[524] | 340 | !------------------------------------------------------------------- |
---|
| 341 | ! --- SET CONSTANTS AND PARAMETERS |
---|
| 342 | !------------------------------------------------------------------- |
---|
[1836] | 343 | ! print *, '-> cv_driver' !jyg |
---|
[524] | 344 | c -- set simulation flags: |
---|
| 345 | c (common cvflag) |
---|
| 346 | |
---|
[1849] | 347 | CALL cv_flag(0) |
---|
[524] | 348 | |
---|
| 349 | c -- set thermodynamical constants: |
---|
| 350 | c (common cvthermo) |
---|
| 351 | |
---|
| 352 | CALL cv_thermo(iflag_con) |
---|
| 353 | |
---|
| 354 | c -- set convect parameters |
---|
| 355 | c |
---|
| 356 | c includes microphysical parameters and parameters that |
---|
| 357 | c control the rate of approach to quasi-equilibrium) |
---|
| 358 | c (common cvparam) |
---|
| 359 | |
---|
[879] | 360 | |
---|
| 361 | if (iflag_con.eq.30) then |
---|
| 362 | CALL cv30_param(nd,delt) |
---|
[524] | 363 | endif |
---|
| 364 | |
---|
| 365 | if (iflag_con.eq.4) then |
---|
| 366 | CALL cv_param(nd) |
---|
| 367 | endif |
---|
| 368 | |
---|
| 369 | !--------------------------------------------------------------------- |
---|
| 370 | ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS |
---|
| 371 | !--------------------------------------------------------------------- |
---|
| 372 | |
---|
[1742] | 373 | ft1(:,:)=0.0 |
---|
| 374 | fq1(:,:)=0.0 |
---|
| 375 | fu1(:,:)=0.0 |
---|
| 376 | fv1(:,:)=0.0 |
---|
| 377 | tvp1(:,:)=0.0 |
---|
| 378 | tp1(:,:)=0.0 |
---|
| 379 | clw1(:,:)=0.0 |
---|
[559] | 380 | cym |
---|
[1742] | 381 | clw(:,:)=0.0 |
---|
| 382 | gz1(:,:) = 0. |
---|
| 383 | VPrecip1(:,:) = 0. |
---|
| 384 | Ma1(:,:)=0.0 |
---|
| 385 | upwd1(:,:)=0.0 |
---|
| 386 | dnwd1(:,:)=0.0 |
---|
| 387 | dnwd01(:,:)=0.0 |
---|
| 388 | qcondc1(:,:)=0.0 |
---|
[524] | 389 | |
---|
[1742] | 390 | ftra1(:,:,:)=0.0 |
---|
[524] | 391 | |
---|
[1742] | 392 | elij1(:,:,:) = 0.0 |
---|
| 393 | sij1(:,:,:) = 0.0 |
---|
[524] | 394 | |
---|
[1742] | 395 | precip1(:)=0.0 |
---|
| 396 | iflag1(:)=0 |
---|
| 397 | wd1(:)=0.0 |
---|
| 398 | cape1(:)=0.0 |
---|
| 399 | |
---|
[879] | 400 | if (iflag_con.eq.30) then |
---|
[524] | 401 | do il=1,len |
---|
| 402 | sig1(il,nd)=sig1(il,nd)+1. |
---|
| 403 | sig1(il,nd)=amin1(sig1(il,nd),12.1) |
---|
| 404 | enddo |
---|
| 405 | endif |
---|
| 406 | |
---|
[1774] | 407 | ! RomP >>> |
---|
| 408 | wdtrainA1(:,:) =0. |
---|
| 409 | wdtrainM1(:,:) =0. |
---|
| 410 | da1(:,:) =0. |
---|
| 411 | phi1(:,:,:) =0. |
---|
| 412 | epmlmMm1(:,:,:) =0. |
---|
| 413 | eplaMm1(:,:) =0. |
---|
| 414 | mp1(:,:) =0. |
---|
| 415 | evap1(:,:) =0. |
---|
| 416 | ep1(:,:) =0. |
---|
| 417 | sij1(:,:,:) =0. |
---|
| 418 | elij1(:,:,:) =0. |
---|
| 419 | phi21(:,:,:) =0. |
---|
| 420 | d1a1(:,:) =0. |
---|
| 421 | dam1(:,:) =0. |
---|
| 422 | ! RomP <<< |
---|
| 423 | |
---|
[524] | 424 | !-------------------------------------------------------------------- |
---|
| 425 | ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY |
---|
| 426 | !-------------------------------------------------------------------- |
---|
| 427 | |
---|
[879] | 428 | if (iflag_con.eq.30) then |
---|
| 429 | |
---|
[965] | 430 | ! print*,'Emanuel version 30 ' |
---|
[879] | 431 | CALL cv30_prelim(len,nd,ndp1,t1,q1,p1,ph1 ! nd->na |
---|
[524] | 432 | o ,lv1,cpn1,tv1,gz1,h1,hm1,th1) |
---|
| 433 | endif |
---|
| 434 | |
---|
| 435 | if (iflag_con.eq.4) then |
---|
| 436 | CALL cv_prelim(len,nd,ndp1,t1,q1,p1,ph1 |
---|
| 437 | o ,lv1,cpn1,tv1,gz1,h1,hm1) |
---|
| 438 | endif |
---|
| 439 | |
---|
| 440 | !-------------------------------------------------------------------- |
---|
| 441 | ! --- CONVECTIVE FEED |
---|
| 442 | !-------------------------------------------------------------------- |
---|
| 443 | |
---|
[879] | 444 | if (iflag_con.eq.30) then |
---|
| 445 | CALL cv30_feed(len,nd,t1,q1,qs1,p1,ph1,hm1,gz1 ! nd->na |
---|
[524] | 446 | o ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1) |
---|
| 447 | endif |
---|
| 448 | |
---|
| 449 | if (iflag_con.eq.4) then |
---|
| 450 | CALL cv_feed(len,nd,t1,q1,qs1,p1,hm1,gz1 |
---|
| 451 | o ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1) |
---|
| 452 | endif |
---|
| 453 | |
---|
| 454 | !-------------------------------------------------------------------- |
---|
| 455 | ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part |
---|
| 456 | ! (up through ICB for convect4, up through ICB+1 for convect3) |
---|
| 457 | ! Calculates the lifted parcel virtual temperature at nk, the |
---|
| 458 | ! actual temperature, and the adiabatic liquid water content. |
---|
| 459 | !-------------------------------------------------------------------- |
---|
| 460 | |
---|
[879] | 461 | if (iflag_con.eq.30) then |
---|
| 462 | CALL cv30_undilute1(len,nd,t1,q1,qs1,gz1,plcl1,p1,nk1,icb1 ! nd->na |
---|
[524] | 463 | o ,tp1,tvp1,clw1,icbs1) |
---|
| 464 | endif |
---|
| 465 | |
---|
| 466 | if (iflag_con.eq.4) then |
---|
| 467 | CALL cv_undilute1(len,nd,t1,q1,qs1,gz1,p1,nk1,icb1,icbmax |
---|
| 468 | : ,tp1,tvp1,clw1) |
---|
| 469 | endif |
---|
| 470 | |
---|
| 471 | !------------------------------------------------------------------- |
---|
| 472 | ! --- TRIGGERING |
---|
| 473 | !------------------------------------------------------------------- |
---|
| 474 | |
---|
[879] | 475 | if (iflag_con.eq.30) then |
---|
| 476 | CALL cv30_trigger(len,nd,icb1,plcl1,p1,th1,tv1,tvp1 ! nd->na |
---|
[524] | 477 | o ,pbase1,buoybase1,iflag1,sig1,w01) |
---|
| 478 | endif |
---|
| 479 | |
---|
| 480 | if (iflag_con.eq.4) then |
---|
| 481 | CALL cv_trigger(len,nd,icb1,cbmf1,tv1,tvp1,iflag1) |
---|
| 482 | endif |
---|
| 483 | |
---|
| 484 | !===================================================================== |
---|
| 485 | ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY |
---|
| 486 | !===================================================================== |
---|
| 487 | |
---|
| 488 | ncum=0 |
---|
| 489 | do 400 i=1,len |
---|
| 490 | if(iflag1(i).eq.0)then |
---|
| 491 | ncum=ncum+1 |
---|
| 492 | idcum(ncum)=i |
---|
| 493 | endif |
---|
| 494 | 400 continue |
---|
| 495 | |
---|
[1836] | 496 | ! print*,'cv_driver : klon, ncum = ',len,ncum |
---|
[524] | 497 | |
---|
| 498 | IF (ncum.gt.0) THEN |
---|
| 499 | |
---|
| 500 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
| 501 | ! --- COMPRESS THE FIELDS |
---|
| 502 | ! (-> vectorization over convective gridpoints) |
---|
| 503 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
| 504 | |
---|
[879] | 505 | if (iflag_con.eq.30) then |
---|
| 506 | CALL cv30_compress( len,nloc,ncum,nd,ntra |
---|
[524] | 507 | : ,iflag1,nk1,icb1,icbs1 |
---|
| 508 | : ,plcl1,tnk1,qnk1,gznk1,pbase1,buoybase1 |
---|
| 509 | : ,t1,q1,qs1,u1,v1,gz1,th1 |
---|
| 510 | : ,tra1 |
---|
| 511 | : ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1 |
---|
| 512 | : ,sig1,w01 |
---|
| 513 | o ,iflag,nk,icb,icbs |
---|
| 514 | o ,plcl,tnk,qnk,gznk,pbase,buoybase |
---|
| 515 | o ,t,q,qs,u,v,gz,th |
---|
| 516 | o ,tra |
---|
| 517 | o ,h,lv,cpn,p,ph,tv,tp,tvp,clw |
---|
| 518 | o ,sig,w0 ) |
---|
| 519 | endif |
---|
| 520 | |
---|
| 521 | if (iflag_con.eq.4) then |
---|
| 522 | CALL cv_compress( len,nloc,ncum,nd |
---|
| 523 | : ,iflag1,nk1,icb1 |
---|
| 524 | : ,cbmf1,plcl1,tnk1,qnk1,gznk1 |
---|
| 525 | : ,t1,q1,qs1,u1,v1,gz1 |
---|
| 526 | : ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1 |
---|
| 527 | o ,iflag,nk,icb |
---|
| 528 | o ,cbmf,plcl,tnk,qnk,gznk |
---|
| 529 | o ,t,q,qs,u,v,gz,h,lv,cpn,p,ph,tv,tp,tvp,clw |
---|
| 530 | o ,dph ) |
---|
| 531 | endif |
---|
| 532 | |
---|
| 533 | !------------------------------------------------------------------- |
---|
| 534 | ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part : |
---|
| 535 | ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
---|
| 536 | ! --- & |
---|
| 537 | ! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE |
---|
| 538 | ! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD |
---|
| 539 | ! --- & |
---|
| 540 | ! --- FIND THE LEVEL OF NEUTRAL BUOYANCY |
---|
| 541 | !------------------------------------------------------------------- |
---|
| 542 | |
---|
[879] | 543 | if (iflag_con.eq.30) then |
---|
| 544 | CALL cv30_undilute2(nloc,ncum,nd,icb,icbs,nk !na->nd |
---|
[524] | 545 | : ,tnk,qnk,gznk,t,q,qs,gz |
---|
| 546 | : ,p,h,tv,lv,pbase,buoybase,plcl |
---|
| 547 | o ,inb,tp,tvp,clw,hp,ep,sigp,buoy) |
---|
| 548 | endif |
---|
| 549 | |
---|
| 550 | if (iflag_con.eq.4) then |
---|
| 551 | CALL cv_undilute2(nloc,ncum,nd,icb,nk |
---|
| 552 | : ,tnk,qnk,gznk,t,q,qs,gz |
---|
| 553 | : ,p,dph,h,tv,lv |
---|
| 554 | o ,inb,inbis,tp,tvp,clw,hp,ep,sigp,frac) |
---|
| 555 | endif |
---|
| 556 | |
---|
| 557 | !------------------------------------------------------------------- |
---|
| 558 | ! --- CLOSURE |
---|
| 559 | !------------------------------------------------------------------- |
---|
| 560 | |
---|
[879] | 561 | if (iflag_con.eq.30) then |
---|
| 562 | CALL cv30_closure(nloc,ncum,nd,icb,inb ! na->nd |
---|
[524] | 563 | : ,pbase,p,ph,tv,buoy |
---|
| 564 | o ,sig,w0,cape,m) |
---|
| 565 | endif |
---|
| 566 | |
---|
| 567 | if (iflag_con.eq.4) then |
---|
| 568 | CALL cv_closure(nloc,ncum,nd,nk,icb |
---|
| 569 | : ,tv,tvp,p,ph,dph,plcl,cpn |
---|
| 570 | o ,iflag,cbmf) |
---|
| 571 | endif |
---|
| 572 | |
---|
| 573 | !------------------------------------------------------------------- |
---|
| 574 | ! --- MIXING |
---|
| 575 | !------------------------------------------------------------------- |
---|
| 576 | |
---|
[879] | 577 | if (iflag_con.eq.30) then |
---|
| 578 | CALL cv30_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb ! na->nd |
---|
[524] | 579 | : ,ph,t,q,qs,u,v,tra,h,lv,qnk |
---|
| 580 | : ,hp,tv,tvp,ep,clw,m,sig |
---|
| 581 | o ,ment,qent,uent,vent,sij,elij,ments,qents,traent) |
---|
| 582 | endif |
---|
| 583 | |
---|
| 584 | if (iflag_con.eq.4) then |
---|
| 585 | CALL cv_mixing(nloc,ncum,nd,icb,nk,inb,inbis |
---|
| 586 | : ,ph,t,q,qs,u,v,h,lv,qnk |
---|
| 587 | : ,hp,tv,tvp,ep,clw,cbmf |
---|
| 588 | o ,m,ment,qent,uent,vent,nent,sij,elij) |
---|
| 589 | endif |
---|
| 590 | |
---|
| 591 | !------------------------------------------------------------------- |
---|
| 592 | ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS |
---|
| 593 | !------------------------------------------------------------------- |
---|
| 594 | |
---|
[879] | 595 | if (iflag_con.eq.30) then |
---|
[1742] | 596 | !RomP >>> |
---|
| 597 | CALL cv30_unsat(nloc,ncum,nd,nd,ntra,icb,inb ! na->nd |
---|
[524] | 598 | : ,t,q,qs,gz,u,v,tra,p,ph |
---|
| 599 | : ,th,tv,lv,cpn,ep,sigp,clw |
---|
| 600 | : ,m,ment,elij,delt,plcl |
---|
[1742] | 601 | : ,mp,qp,up,vp,trap,wt,water,evap,b |
---|
| 602 | o ,wdtrainA,wdtrainM) |
---|
| 603 | !RomP <<< |
---|
[524] | 604 | endif |
---|
| 605 | |
---|
| 606 | if (iflag_con.eq.4) then |
---|
| 607 | CALL cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph |
---|
| 608 | : ,h,lv,ep,sigp,clw,m,ment,elij |
---|
| 609 | o ,iflag,mp,qp,up,vp,wt,water,evap) |
---|
| 610 | endif |
---|
| 611 | |
---|
| 612 | !------------------------------------------------------------------- |
---|
| 613 | ! --- YIELD |
---|
| 614 | ! (tendencies, precipitation, variables of interface with other |
---|
| 615 | ! processes, etc) |
---|
| 616 | !------------------------------------------------------------------- |
---|
| 617 | |
---|
[879] | 618 | if (iflag_con.eq.30) then |
---|
| 619 | CALL cv30_yield(nloc,ncum,nd,nd,ntra ! na->nd |
---|
[524] | 620 | : ,icb,inb,delt |
---|
| 621 | : ,t,q,u,v,tra,gz,p,ph,h,hp,lv,cpn,th |
---|
| 622 | : ,ep,clw,m,tp,mp,qp,up,vp,trap |
---|
| 623 | : ,wt,water,evap,b |
---|
| 624 | : ,ment,qent,uent,vent,nent,elij,traent,sig |
---|
| 625 | : ,tv,tvp |
---|
[619] | 626 | o ,iflag,precip,VPrecip,ft,fq,fu,fv,ftra |
---|
[524] | 627 | o ,upwd,dnwd,dnwd0,ma,mike,tls,tps,qcondc,wd) |
---|
| 628 | endif |
---|
| 629 | |
---|
| 630 | if (iflag_con.eq.4) then |
---|
| 631 | CALL cv_yield(nloc,ncum,nd,nk,icb,inb,delt |
---|
| 632 | : ,t,q,u,v,gz,p,ph,h,hp,lv,cpn |
---|
| 633 | : ,ep,clw,frac,m,mp,qp,up,vp |
---|
| 634 | : ,wt,water,evap |
---|
| 635 | : ,ment,qent,uent,vent,nent,elij |
---|
| 636 | : ,tv,tvp |
---|
| 637 | o ,iflag,wd,qprime,tprime |
---|
| 638 | o ,precip,cbmf,ft,fq,fu,fv,Ma,qcondc) |
---|
| 639 | endif |
---|
| 640 | |
---|
| 641 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
[619] | 642 | ! --- passive tracers |
---|
[524] | 643 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
| 644 | |
---|
[879] | 645 | if (iflag_con.eq.30) then |
---|
[1742] | 646 | !RomP >>> |
---|
[879] | 647 | CALL cv30_tracer(nloc,len,ncum,nd,nd, |
---|
[1742] | 648 | : ment,sij,da,phi,phi2,d1a,dam, |
---|
| 649 | : ep,VPrecip,elij,clw,epmlmMm,eplaMm, |
---|
| 650 | : icb,inb) |
---|
| 651 | !RomP <<< |
---|
[619] | 652 | endif |
---|
[524] | 653 | |
---|
[619] | 654 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
| 655 | ! --- UNCOMPRESS THE FIELDS |
---|
| 656 | !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ |
---|
| 657 | c set iflag1 =42 for non convective points |
---|
| 658 | do i=1,len |
---|
| 659 | iflag1(i)=42 |
---|
| 660 | end do |
---|
| 661 | c |
---|
[879] | 662 | if (iflag_con.eq.30) then |
---|
| 663 | CALL cv30_uncompress(nloc,len,ncum,nd,ntra,idcum |
---|
[524] | 664 | : ,iflag |
---|
[1742] | 665 | : ,precip,VPrecip,evap,ep,sig,w0 !RomP |
---|
[524] | 666 | : ,ft,fq,fu,fv,ftra |
---|
[1742] | 667 | : ,inb |
---|
[524] | 668 | : ,Ma,upwd,dnwd,dnwd0,qcondc,wd,cape |
---|
[1742] | 669 | : ,da,phi,mp,phi2,d1a,dam,sij !RomP |
---|
| 670 | : ,elij,clw,epmlmMm,eplaMm !RomP |
---|
| 671 | : ,wdtrainA,wdtrainM !RomP |
---|
[524] | 672 | o ,iflag1 |
---|
[1742] | 673 | o ,precip1,VPrecip1,evap1,ep1,sig1,w01 !RomP |
---|
[524] | 674 | o ,ft1,fq1,fu1,fv1,ftra1 |
---|
[619] | 675 | o ,inb1 |
---|
| 676 | o ,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1 |
---|
[1742] | 677 | o ,da1,phi1,mp1,phi21,d1a1,dam1,sij1 !RomP |
---|
| 678 | o ,elij1,clw1,epmlmMm1,eplaMm1 !RomP |
---|
| 679 | o ,wdtrainA1,wdtrainM1) !RomP |
---|
[524] | 680 | endif |
---|
| 681 | |
---|
| 682 | if (iflag_con.eq.4) then |
---|
| 683 | CALL cv_uncompress(nloc,len,ncum,nd,idcum |
---|
| 684 | : ,iflag |
---|
| 685 | : ,precip,cbmf |
---|
| 686 | : ,ft,fq,fu,fv |
---|
| 687 | : ,Ma,qcondc |
---|
| 688 | o ,iflag1 |
---|
| 689 | o ,precip1,cbmf1 |
---|
| 690 | o ,ft1,fq1,fu1,fv1 |
---|
| 691 | o ,Ma1,qcondc1 ) |
---|
| 692 | endif |
---|
| 693 | |
---|
| 694 | ENDIF ! ncum>0 |
---|
| 695 | |
---|
| 696 | 9999 continue |
---|
| 697 | |
---|
[1836] | 698 | ! print *, 'fin cv_driver ->' !jyg |
---|
[524] | 699 | return |
---|
| 700 | end |
---|
| 701 | |
---|
| 702 | !================================================================== |
---|
[1849] | 703 | SUBROUTINE cv_flag(iflag_ice_thermo) |
---|
[524] | 704 | implicit none |
---|
| 705 | |
---|
[1849] | 706 | c Argument : iflag_ice_thermo : ice thermodynamics is taken into account if |
---|
| 707 | c iflag_ice_thermo >=1 |
---|
| 708 | INTEGER iflag_ice_thermo |
---|
| 709 | |
---|
[524] | 710 | #include "cvflag.h" |
---|
| 711 | |
---|
| 712 | c -- si .TRUE., on rend la gravite plus explicite et eventuellement |
---|
| 713 | c differente de 10.0 dans convect3: |
---|
[619] | 714 | cvflag_grav = .TRUE. |
---|
[1849] | 715 | cvflag_ice = iflag_ice_thermo .GE. 1 |
---|
[524] | 716 | |
---|
| 717 | return |
---|
| 718 | end |
---|
| 719 | |
---|
| 720 | !================================================================== |
---|
| 721 | SUBROUTINE cv_thermo(iflag_con) |
---|
| 722 | implicit none |
---|
| 723 | |
---|
| 724 | c------------------------------------------------------------- |
---|
| 725 | c Set thermodynamical constants for convectL |
---|
| 726 | c------------------------------------------------------------- |
---|
| 727 | |
---|
| 728 | #include "YOMCST.h" |
---|
| 729 | #include "cvthermo.h" |
---|
| 730 | |
---|
| 731 | integer iflag_con |
---|
| 732 | |
---|
| 733 | |
---|
| 734 | c original set from convect: |
---|
| 735 | if (iflag_con.eq.4) then |
---|
| 736 | cpd=1005.7 |
---|
| 737 | cpv=1870.0 |
---|
| 738 | cl=4190.0 |
---|
| 739 | rrv=461.5 |
---|
| 740 | rrd=287.04 |
---|
| 741 | lv0=2.501E6 |
---|
| 742 | g=9.8 |
---|
| 743 | t0=273.15 |
---|
| 744 | grav=g |
---|
[879] | 745 | else |
---|
[524] | 746 | |
---|
| 747 | c constants consistent with LMDZ: |
---|
| 748 | cpd = RCPD |
---|
| 749 | cpv = RCPV |
---|
| 750 | cl = RCW |
---|
[1849] | 751 | ci = RCS |
---|
[524] | 752 | rrv = RV |
---|
| 753 | rrd = RD |
---|
| 754 | lv0 = RLVTT |
---|
[1849] | 755 | lf0 = RLSTT-RLVTT |
---|
[524] | 756 | g = RG ! not used in convect3 |
---|
| 757 | c ori t0 = RTT |
---|
| 758 | t0 = 273.15 ! convect3 (RTT=273.16) |
---|
[619] | 759 | c maf grav= 10. ! implicitely or explicitely used in convect3 |
---|
| 760 | grav= g ! implicitely or explicitely used in convect3 |
---|
[524] | 761 | endif |
---|
| 762 | |
---|
| 763 | rowl=1000.0 !(a quelle variable de YOMCST cela correspond-il?) |
---|
| 764 | |
---|
| 765 | clmcpv=cl-cpv |
---|
| 766 | clmcpd=cl-cpd |
---|
[1849] | 767 | clmci=cl-ci |
---|
[524] | 768 | cpdmcp=cpd-cpv |
---|
| 769 | cpvmcpd=cpv-cpd |
---|
| 770 | cpvmcl=cl-cpv ! for convect3 |
---|
| 771 | eps=rrd/rrv |
---|
| 772 | epsi=1.0/eps |
---|
| 773 | epsim1=epsi-1.0 |
---|
| 774 | c ginv=1.0/g |
---|
| 775 | ginv=1.0/grav |
---|
| 776 | hrd=0.5*rrd |
---|
| 777 | |
---|
| 778 | return |
---|
[1849] | 779 | end |
---|