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