[4590] | 1 | MODULE lmdz_thermcell_main |
---|
[1403] | 2 | ! $Id: lmdz_thermcell_main.F90 4682 2023-09-08 10:08:46Z fhourdin $ |
---|
[878] | 3 | ! |
---|
[4678] | 4 | ! A REGARDER !!!!!!!!!!!!!!!!! |
---|
| 5 | ! ATTENTION : zpspsk est inout et out mais c'est pas forcement pour de bonnes raisons (FH, 2023) |
---|
[4590] | 6 | CONTAINS |
---|
| 7 | |
---|
[4089] | 8 | subroutine thermcell_main(itap,ngrid,nlay,ptimestep & |
---|
[878] | 9 | & ,pplay,pplev,pphi,debut & |
---|
[4678] | 10 | & ,pu,pv,pt,p_o & |
---|
[878] | 11 | & ,pduadj,pdvadj,pdtadj,pdoadj & |
---|
[1026] | 12 | & ,fm0,entr0,detr0,zqta,zqla,lmax & |
---|
[878] | 13 | & ,ratqscth,ratqsdiff,zqsatth & |
---|
[1403] | 14 | & ,zmax0, f0,zw2,fraca,ztv & |
---|
[4089] | 15 | & ,zpspsk,ztla,zthl,ztva & |
---|
| 16 | & ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax & |
---|
| 17 | #ifdef ISO |
---|
| 18 | & ,xtpo,xtpdoadj & |
---|
| 19 | #endif |
---|
| 20 | & ) |
---|
[878] | 21 | |
---|
[4089] | 22 | |
---|
[4590] | 23 | USE lmdz_thermcell_ini, ONLY: thermcell_ini,dqimpl,dvdq,prt_level,lunout,prt_level |
---|
| 24 | USE lmdz_thermcell_ini, ONLY: iflag_thermals_closure,iflag_thermals_ed,tau_thermals,r_aspect_thermals |
---|
| 25 | USE lmdz_thermcell_ini, ONLY: iflag_thermals_down,fact_thermals_down |
---|
| 26 | USE lmdz_thermcell_ini, ONLY: RD,RG |
---|
[4089] | 27 | |
---|
[4590] | 28 | USE lmdz_thermcell_down, ONLY: thermcell_updown_dq |
---|
| 29 | USE lmdz_thermcell_closure, ONLY: thermcell_closure |
---|
| 30 | USE lmdz_thermcell_dq, ONLY: thermcell_dq |
---|
| 31 | USE lmdz_thermcell_dry, ONLY: thermcell_dry |
---|
| 32 | USE lmdz_thermcell_dv2, ONLY: thermcell_dv2 |
---|
| 33 | USE lmdz_thermcell_env, ONLY: thermcell_env |
---|
| 34 | USE lmdz_thermcell_flux2, ONLY: thermcell_flux2 |
---|
| 35 | USE lmdz_thermcell_height, ONLY: thermcell_height |
---|
| 36 | USE lmdz_thermcell_plume, ONLY: thermcell_plume |
---|
| 37 | USE lmdz_thermcell_plume_6A, ONLY: thermcell_plume_6A,thermcell_plume_5B |
---|
| 38 | |
---|
[4089] | 39 | #ifdef ISO |
---|
[4143] | 40 | USE infotrac_phy, ONLY : ntiso |
---|
[4089] | 41 | #ifdef ISOVERIF |
---|
| 42 | USE isotopes_mod, ONLY : iso_eau,iso_HDO |
---|
| 43 | USE isotopes_verif_mod, ONLY: iso_verif_egalite, & |
---|
| 44 | iso_verif_aberrant_encadre |
---|
| 45 | #endif |
---|
| 46 | #endif |
---|
| 47 | |
---|
| 48 | |
---|
[878] | 49 | IMPLICIT NONE |
---|
| 50 | |
---|
| 51 | !======================================================================= |
---|
| 52 | ! Auteurs: Frederic Hourdin, Catherine Rio, Anne Mathieu |
---|
| 53 | ! Version du 09.02.07 |
---|
| 54 | ! Calcul du transport vertical dans la couche limite en presence |
---|
| 55 | ! de "thermiques" explicitement representes avec processus nuageux |
---|
| 56 | ! |
---|
[1403] | 57 | ! Reecriture a partir d'un listing papier a Habas, le 14/02/00 |
---|
[878] | 58 | ! |
---|
[1403] | 59 | ! le thermique est suppose homogene et dissipe par melange avec |
---|
| 60 | ! son environnement. la longueur l_mix controle l'efficacite du |
---|
| 61 | ! melange |
---|
[878] | 62 | ! |
---|
[1403] | 63 | ! Le calcul du transport des differentes especes se fait en prenant |
---|
[878] | 64 | ! en compte: |
---|
| 65 | ! 1. un flux de masse montant |
---|
| 66 | ! 2. un flux de masse descendant |
---|
| 67 | ! 3. un entrainement |
---|
| 68 | ! 4. un detrainement |
---|
| 69 | ! |
---|
[1738] | 70 | ! Modif 2013/01/04 (FH hourdin@lmd.jussieu.fr) |
---|
| 71 | ! Introduction of an implicit computation of vertical advection in |
---|
| 72 | ! the environment of thermal plumes in thermcell_dq |
---|
| 73 | ! impl = 0 : explicit, 1 : implicit, -1 : old version |
---|
| 74 | ! controled by iflag_thermals = |
---|
| 75 | ! 15, 16 run with impl=-1 : numerical convergence with NPv3 |
---|
| 76 | ! 17, 18 run with impl=1 : more stable |
---|
| 77 | ! 15 and 17 correspond to the activation of the stratocumulus "bidouille" |
---|
| 78 | ! |
---|
[4133] | 79 | ! Using |
---|
| 80 | ! abort_physic |
---|
| 81 | ! iso_verif_aberrant_encadre |
---|
| 82 | ! iso_verif_egalite |
---|
| 83 | ! test_ltherm |
---|
| 84 | ! thermcell_closure |
---|
| 85 | ! thermcell_dq |
---|
| 86 | ! thermcell_dry |
---|
| 87 | ! thermcell_dv2 |
---|
| 88 | ! thermcell_env |
---|
| 89 | ! thermcell_flux2 |
---|
| 90 | ! thermcell_height |
---|
| 91 | ! thermcell_plume |
---|
| 92 | ! thermcell_plume_5B |
---|
| 93 | ! thermcell_plume_6A |
---|
| 94 | ! |
---|
[878] | 95 | !======================================================================= |
---|
| 96 | |
---|
[1738] | 97 | |
---|
[878] | 98 | !----------------------------------------------------------------------- |
---|
| 99 | ! declarations: |
---|
| 100 | ! ------------- |
---|
| 101 | |
---|
| 102 | |
---|
| 103 | ! arguments: |
---|
| 104 | ! ---------- |
---|
[4089] | 105 | integer, intent(in) :: itap,ngrid,nlay |
---|
| 106 | real, intent(in) :: ptimestep |
---|
[4590] | 107 | real, intent(in), dimension(ngrid,nlay) :: pt,pu,pv,pplay,pphi |
---|
[4678] | 108 | ! ATTENTION : zpspsk est inout et out mais c'est pas forcement pour de bonnes raisons (FH, 2023) |
---|
| 109 | real, intent(in), dimension(ngrid,nlay) :: p_o |
---|
[4590] | 110 | real, intent(out), dimension(ngrid,nlay) :: zpspsk |
---|
[4089] | 111 | real, intent(in), dimension(ngrid,nlay+1) :: pplev |
---|
[4133] | 112 | integer, intent(out), dimension(ngrid) :: lmax |
---|
[4089] | 113 | real, intent(out), dimension(ngrid,nlay) :: pdtadj,pduadj,pdvadj,pdoadj,entr0,detr0 |
---|
| 114 | real, intent(out), dimension(ngrid,nlay) :: ztla,zqla,zqta,zqsatth,zthl |
---|
| 115 | real, intent(out), dimension(ngrid,nlay+1) :: fm0,zw2,fraca |
---|
[4133] | 116 | real, intent(inout), dimension(ngrid) :: zmax0,f0 |
---|
[4089] | 117 | real, intent(out), dimension(ngrid,nlay) :: ztva,ztv |
---|
| 118 | logical, intent(in) :: debut |
---|
[4133] | 119 | real,intent(out), dimension(ngrid,nlay) :: ratqscth,ratqsdiff |
---|
[878] | 120 | |
---|
[4089] | 121 | real, intent(out), dimension(ngrid) :: pcon |
---|
| 122 | real, intent(out), dimension(ngrid,nlay) :: rhobarz,wth3 |
---|
| 123 | real, intent(out), dimension(ngrid) :: wmax_sec |
---|
| 124 | integer,intent(out), dimension(ngrid) :: lalim |
---|
| 125 | real, intent(out), dimension(ngrid,nlay+1) :: fm |
---|
| 126 | real, intent(out), dimension(ngrid,nlay) :: alim_star |
---|
| 127 | real, intent(out), dimension(ngrid) :: zmax |
---|
[972] | 128 | |
---|
[878] | 129 | ! local: |
---|
| 130 | ! ------ |
---|
| 131 | |
---|
[1738] | 132 | |
---|
[883] | 133 | integer,save :: igout=1 |
---|
[987] | 134 | !$OMP THREADPRIVATE(igout) |
---|
[938] | 135 | integer,save :: lunout1=6 |
---|
[987] | 136 | !$OMP THREADPRIVATE(lunout1) |
---|
[883] | 137 | integer,save :: lev_out=10 |
---|
[987] | 138 | !$OMP THREADPRIVATE(lev_out) |
---|
[878] | 139 | |
---|
[4089] | 140 | real lambda, zf,zf2,var,vardiff,CHI |
---|
| 141 | integer ig,k,l,ierr,ll |
---|
[878] | 142 | logical sorties |
---|
[4089] | 143 | real, dimension(ngrid) :: linter,zmix, zmax_sec |
---|
[4133] | 144 | integer,dimension(ngrid) :: lmin,lmix,lmix_bis,nivcon |
---|
[4089] | 145 | real, dimension(ngrid,nlay) :: ztva_est |
---|
[4678] | 146 | real, dimension(ngrid,nlay) :: deltaz,zlay,zh,zdthladj,zu,zv,z_o,zl,zva,zua,z_oa |
---|
[4089] | 147 | real, dimension(ngrid,nlay) :: zta,zha,q2,wq,wthl,wthv,thetath2,wth2 |
---|
[4133] | 148 | real, dimension(ngrid,nlay) :: rho,masse |
---|
[4089] | 149 | real, dimension(ngrid,nlay+1) :: zw_est,zlev |
---|
| 150 | real, dimension(ngrid) :: wmax,wmax_tmp |
---|
| 151 | real, dimension(ngrid,nlay+1) :: f_star |
---|
| 152 | real, dimension(ngrid,nlay) :: entr,detr,entr_star,detr_star,alim_star_clos |
---|
| 153 | real, dimension(ngrid,nlay) :: zqsat,csc |
---|
| 154 | real, dimension(ngrid) :: zcon,zcon2,alim_star_tot,f |
---|
[4413] | 155 | real, dimension(ngrid,nlay) :: entrdn,detrdn |
---|
[878] | 156 | |
---|
[4089] | 157 | character (len=20) :: modname='thermcell_main' |
---|
| 158 | character (len=80) :: abort_message |
---|
[878] | 159 | |
---|
| 160 | |
---|
[4089] | 161 | #ifdef ISO |
---|
[4143] | 162 | REAL xtpo(ntiso,ngrid,nlay),xtpdoadj(ntiso,ngrid,nlay) |
---|
| 163 | REAL xtzo(ntiso,ngrid,nlay) |
---|
[4089] | 164 | REAL xtpdoadj_tmp(ngrid,nlay) |
---|
| 165 | REAL xtpo_tmp(ngrid,nlay) |
---|
| 166 | REAL xtzo_tmp(ngrid,nlay) |
---|
| 167 | integer ixt |
---|
| 168 | #endif |
---|
[878] | 169 | |
---|
| 170 | ! |
---|
| 171 | |
---|
| 172 | !----------------------------------------------------------------------- |
---|
| 173 | ! initialisation: |
---|
| 174 | ! --------------- |
---|
| 175 | ! |
---|
[1943] | 176 | fm=0. ; entr=0. ; detr=0. |
---|
[972] | 177 | |
---|
[938] | 178 | if (prt_level.ge.1) print*,'thermcell_main V4' |
---|
[878] | 179 | |
---|
| 180 | sorties=.true. |
---|
[4089] | 181 | IF(ngrid.NE.ngrid) THEN |
---|
[878] | 182 | PRINT* |
---|
| 183 | PRINT*,'STOP dans convadj' |
---|
| 184 | PRINT*,'ngrid =',ngrid |
---|
[4089] | 185 | PRINT*,'ngrid =',ngrid |
---|
[878] | 186 | ENDIF |
---|
| 187 | ! |
---|
[1403] | 188 | ! write(lunout,*)'WARNING thermcell_main f0=max(f0,1.e-2)' |
---|
[4089] | 189 | do ig=1,ngrid |
---|
[972] | 190 | f0(ig)=max(f0(ig),1.e-2) |
---|
[1403] | 191 | zmax0(ig)=max(zmax0(ig),40.) |
---|
[972] | 192 | !IMmarche pas ?! if (f0(ig)<1.e-2) f0(ig)=1.e-2 |
---|
| 193 | enddo |
---|
[878] | 194 | |
---|
[1494] | 195 | if (prt_level.ge.20) then |
---|
| 196 | do ig=1,ngrid |
---|
| 197 | print*,'th_main ig f0',ig,f0(ig) |
---|
| 198 | enddo |
---|
| 199 | endif |
---|
[878] | 200 | !----------------------------------------------------------------------- |
---|
| 201 | ! Calcul de T,q,ql a partir de Tl et qT dans l environnement |
---|
| 202 | ! -------------------------------------------------------------------- |
---|
| 203 | ! |
---|
[4678] | 204 | CALL thermcell_env(ngrid,nlay,p_o,pt,pu,pv,pplay, & |
---|
| 205 | & pplev,z_o,zh,zl,ztv,zthl,zu,zv,zpspsk,zqsat,lev_out) |
---|
[878] | 206 | |
---|
[938] | 207 | if (prt_level.ge.1) print*,'thermcell_main apres thermcell_env' |
---|
[878] | 208 | |
---|
| 209 | !------------------------------------------------------------------------ |
---|
| 210 | ! -------------------- |
---|
| 211 | ! |
---|
| 212 | ! |
---|
| 213 | ! + + + + + + + + + + + |
---|
| 214 | ! |
---|
| 215 | ! |
---|
| 216 | ! wa, fraca, wd, fracd -------------------- zlev(2), rhobarz |
---|
| 217 | ! wh,wt,wo ... |
---|
| 218 | ! |
---|
[4678] | 219 | ! + + + + + + + + + + + zh,zu,zv,z_o,rho |
---|
[878] | 220 | ! |
---|
| 221 | ! |
---|
| 222 | ! -------------------- zlev(1) |
---|
| 223 | ! \\\\\\\\\\\\\\\\\\\\ |
---|
| 224 | ! |
---|
| 225 | ! |
---|
| 226 | |
---|
| 227 | !----------------------------------------------------------------------- |
---|
| 228 | ! Calcul des altitudes des couches |
---|
| 229 | !----------------------------------------------------------------------- |
---|
| 230 | |
---|
| 231 | do l=2,nlay |
---|
| 232 | zlev(:,l)=0.5*(pphi(:,l)+pphi(:,l-1))/RG |
---|
| 233 | enddo |
---|
[4089] | 234 | zlev(:,1)=0. |
---|
| 235 | zlev(:,nlay+1)=(2.*pphi(:,nlay)-pphi(:,nlay-1))/RG |
---|
[878] | 236 | do l=1,nlay |
---|
| 237 | zlay(:,l)=pphi(:,l)/RG |
---|
| 238 | enddo |
---|
| 239 | do l=1,nlay |
---|
| 240 | deltaz(:,l)=zlev(:,l+1)-zlev(:,l) |
---|
| 241 | enddo |
---|
| 242 | |
---|
| 243 | !----------------------------------------------------------------------- |
---|
[4089] | 244 | ! Calcul des densites et masses |
---|
[878] | 245 | !----------------------------------------------------------------------- |
---|
| 246 | |
---|
[4089] | 247 | rho(:,:)=pplay(:,:)/(zpspsk(:,:)*RD*ztv(:,:)) |
---|
| 248 | if (prt_level.ge.10) write(lunout,*) 'WARNING thermcell_main rhobarz(:,1)=rho(:,1)' |
---|
[972] | 249 | rhobarz(:,1)=rho(:,1) |
---|
[878] | 250 | do l=2,nlay |
---|
| 251 | rhobarz(:,l)=0.5*(rho(:,l)+rho(:,l-1)) |
---|
| 252 | enddo |
---|
| 253 | do l=1,nlay |
---|
| 254 | masse(:,l)=(pplev(:,l)-pplev(:,l+1))/RG |
---|
| 255 | enddo |
---|
[938] | 256 | if (prt_level.ge.1) print*,'thermcell_main apres initialisation' |
---|
[878] | 257 | |
---|
| 258 | !------------------------------------------------------------------ |
---|
| 259 | ! |
---|
| 260 | ! /|\ |
---|
| 261 | ! -------- | F_k+1 ------- |
---|
| 262 | ! ----> D_k |
---|
| 263 | ! /|\ <---- E_k , A_k |
---|
| 264 | ! -------- | F_k --------- |
---|
| 265 | ! ----> D_k-1 |
---|
| 266 | ! <---- E_k-1 , A_k-1 |
---|
| 267 | ! |
---|
| 268 | ! |
---|
| 269 | ! |
---|
| 270 | ! |
---|
| 271 | ! |
---|
| 272 | ! --------------------------- |
---|
| 273 | ! |
---|
| 274 | ! ----- F_lmax+1=0 ---------- \ |
---|
| 275 | ! lmax (zmax) | |
---|
| 276 | ! --------------------------- | |
---|
| 277 | ! | |
---|
| 278 | ! --------------------------- | |
---|
| 279 | ! | |
---|
| 280 | ! --------------------------- | |
---|
| 281 | ! | |
---|
| 282 | ! --------------------------- | |
---|
| 283 | ! | |
---|
| 284 | ! --------------------------- | |
---|
| 285 | ! | E |
---|
| 286 | ! --------------------------- | D |
---|
| 287 | ! | |
---|
| 288 | ! --------------------------- | |
---|
| 289 | ! | |
---|
| 290 | ! --------------------------- \ | |
---|
| 291 | ! lalim | | |
---|
| 292 | ! --------------------------- | | |
---|
| 293 | ! | | |
---|
| 294 | ! --------------------------- | | |
---|
| 295 | ! | A | |
---|
| 296 | ! --------------------------- | | |
---|
| 297 | ! | | |
---|
| 298 | ! --------------------------- | | |
---|
| 299 | ! lmin (=1 pour le moment) | | |
---|
| 300 | ! ----- F_lmin=0 ------------ / / |
---|
| 301 | ! |
---|
| 302 | ! --------------------------- |
---|
| 303 | ! ////////////////////////// |
---|
| 304 | ! |
---|
| 305 | ! |
---|
| 306 | !============================================================================= |
---|
| 307 | ! Calculs initiaux ne faisant pas intervenir les changements de phase |
---|
| 308 | !============================================================================= |
---|
| 309 | |
---|
| 310 | !------------------------------------------------------------------ |
---|
[1403] | 311 | ! 1. alim_star est le profil vertical de l'alimentation a la base du |
---|
| 312 | ! panache thermique, calcule a partir de la flotabilite de l'air sec |
---|
[878] | 313 | ! 2. lmin et lalim sont les indices inferieurs et superieurs de alim_star |
---|
| 314 | !------------------------------------------------------------------ |
---|
| 315 | ! |
---|
| 316 | entr_star=0. ; detr_star=0. ; alim_star=0. ; alim_star_tot=0. |
---|
[1403] | 317 | lmin=1 |
---|
[878] | 318 | |
---|
| 319 | !----------------------------------------------------------------------------- |
---|
| 320 | ! 3. wmax_sec et zmax_sec sont les vitesses et altitudes maximum d'un |
---|
| 321 | ! panache sec conservatif (e=d=0) alimente selon alim_star |
---|
| 322 | ! Il s'agit d'un calcul de type CAPE |
---|
[1403] | 323 | ! zmax_sec est utilise pour determiner la geometrie du thermique. |
---|
[878] | 324 | !------------------------------------------------------------------------------ |
---|
[1403] | 325 | !--------------------------------------------------------------------------------- |
---|
| 326 | !calcul du melange et des variables dans le thermique |
---|
| 327 | !-------------------------------------------------------------------------------- |
---|
[878] | 328 | ! |
---|
[1403] | 329 | if (prt_level.ge.1) print*,'avant thermcell_plume ',lev_out |
---|
[878] | 330 | |
---|
[3451] | 331 | !===================================================================== |
---|
| 332 | ! Old version of thermcell_plume in thermcell_plume_6A.F90 |
---|
| 333 | ! It includes both thermcell_plume_6A and thermcell_plume_5B corresponding |
---|
| 334 | ! to the 5B and 6A versions used for CMIP5 and CMIP6. |
---|
| 335 | ! The latest was previously named thermcellV1_plume. |
---|
| 336 | ! The new thermcell_plume is a clean version (removing obsolete |
---|
| 337 | ! options) of thermcell_plume_6A. |
---|
| 338 | ! The 3 versions are controled by |
---|
| 339 | ! flag_thermals_ed <= 9 thermcell_plume_6A |
---|
| 340 | ! <= 19 thermcell_plume_5B |
---|
| 341 | ! else thermcell_plume (default 20 for convergence with 6A) |
---|
| 342 | ! Fredho |
---|
| 343 | !===================================================================== |
---|
[878] | 344 | |
---|
[1403] | 345 | if (iflag_thermals_ed<=9) then |
---|
| 346 | ! print*,'THERM NOUVELLE/NOUVELLE Arnaud' |
---|
[4678] | 347 | CALL thermcell_plume_6A(itap,ngrid,nlay,ptimestep,ztv,zthl,p_o,zl,rhobarz,& |
---|
[1403] | 348 | & zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot, & |
---|
| 349 | & lalim,f0,detr_star,entr_star,f_star,csc,ztva, & |
---|
| 350 | & ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter & |
---|
| 351 | & ,lev_out,lunout1,igout) |
---|
[878] | 352 | |
---|
[3451] | 353 | elseif (iflag_thermals_ed<=19) then |
---|
[1403] | 354 | ! print*,'THERM RIO et al 2010, version d Arnaud' |
---|
[4678] | 355 | CALL thermcell_plume_5B(itap,ngrid,nlay,ptimestep,ztv,zthl,p_o,zl,rhobarz,& |
---|
[1403] | 356 | & zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot, & |
---|
| 357 | & lalim,f0,detr_star,entr_star,f_star,csc,ztva, & |
---|
| 358 | & ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter & |
---|
| 359 | & ,lev_out,lunout1,igout) |
---|
[3451] | 360 | else |
---|
[4678] | 361 | CALL thermcell_plume(itap,ngrid,nlay,ptimestep,ztv,zthl,p_o,zl,rhobarz,& |
---|
[3451] | 362 | & zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot, & |
---|
| 363 | & lalim,f0,detr_star,entr_star,f_star,csc,ztva, & |
---|
| 364 | & ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter & |
---|
| 365 | & ,lev_out,lunout1,igout) |
---|
[1403] | 366 | endif |
---|
[878] | 367 | |
---|
[972] | 368 | if (prt_level.ge.1) print*,'apres thermcell_plume ',lev_out |
---|
| 369 | |
---|
[4678] | 370 | call test_ltherm(ngrid,nlay,pplay,lalim,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_plum lalim ') |
---|
| 371 | call test_ltherm(ngrid,nlay,pplay,lmix ,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_plum lmix ') |
---|
[878] | 372 | |
---|
[938] | 373 | if (prt_level.ge.1) print*,'thermcell_main apres thermcell_plume' |
---|
| 374 | if (prt_level.ge.10) then |
---|
[972] | 375 | write(lunout1,*) 'Dans thermcell_main 2' |
---|
| 376 | write(lunout1,*) 'lmin ',lmin(igout) |
---|
| 377 | write(lunout1,*) 'lalim ',lalim(igout) |
---|
| 378 | write(lunout1,*) ' ig l alim_star entr_star detr_star f_star ' |
---|
| 379 | write(lunout1,'(i6,i4,4e15.5)') (igout,l,alim_star(igout,l),entr_star(igout,l),detr_star(igout,l) & |
---|
[878] | 380 | & ,f_star(igout,l+1),l=1,nint(linter(igout))+5) |
---|
| 381 | endif |
---|
| 382 | |
---|
| 383 | !------------------------------------------------------------------------------- |
---|
| 384 | ! Calcul des caracteristiques du thermique:zmax,zmix,wmax |
---|
| 385 | !------------------------------------------------------------------------------- |
---|
| 386 | ! |
---|
| 387 | CALL thermcell_height(ngrid,nlay,lalim,lmin,linter,lmix,zw2, & |
---|
[4094] | 388 | & zlev,lmax,zmax,zmax0,zmix,wmax) |
---|
[1403] | 389 | ! Attention, w2 est transforme en sa racine carree dans cette routine |
---|
[4143] | 390 | ! Le probleme vient du fait que linter et lmix sont souvent egaux a 1. |
---|
[1403] | 391 | wmax_tmp=0. |
---|
| 392 | do l=1,nlay |
---|
| 393 | wmax_tmp(:)=max(wmax_tmp(:),zw2(:,l)) |
---|
| 394 | enddo |
---|
| 395 | ! print*,"ZMAX ",lalim,lmin,linter,lmix,lmax,zmax,zmax0,zmix,wmax |
---|
[878] | 396 | |
---|
| 397 | |
---|
[1403] | 398 | |
---|
[4678] | 399 | call test_ltherm(ngrid,nlay,pplay,lalim,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_heig lalim ') |
---|
| 400 | call test_ltherm(ngrid,nlay,pplay,lmin ,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_heig lmin ') |
---|
| 401 | call test_ltherm(ngrid,nlay,pplay,lmix ,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_heig lmix ') |
---|
| 402 | call test_ltherm(ngrid,nlay,pplay,lmax ,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_heig lmax ') |
---|
[878] | 403 | |
---|
[938] | 404 | if (prt_level.ge.1) print*,'thermcell_main apres thermcell_height' |
---|
[878] | 405 | |
---|
| 406 | !------------------------------------------------------------------------------- |
---|
| 407 | ! Fermeture,determination de f |
---|
| 408 | !------------------------------------------------------------------------------- |
---|
[1026] | 409 | ! |
---|
[1403] | 410 | ! |
---|
| 411 | CALL thermcell_dry(ngrid,nlay,zlev,pphi,ztv,alim_star, & |
---|
[4094] | 412 | & lalim,lmin,zmax_sec,wmax_sec) |
---|
[878] | 413 | |
---|
[1998] | 414 | |
---|
[4678] | 415 | call test_ltherm(ngrid,nlay,pplay,lmin,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_dry lmin ') |
---|
| 416 | call test_ltherm(ngrid,nlay,pplay,lalim,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_dry lalim ') |
---|
[1403] | 417 | |
---|
| 418 | if (prt_level.ge.1) print*,'thermcell_main apres thermcell_dry' |
---|
| 419 | if (prt_level.ge.10) then |
---|
| 420 | write(lunout1,*) 'Dans thermcell_main 1b' |
---|
| 421 | write(lunout1,*) 'lmin ',lmin(igout) |
---|
| 422 | write(lunout1,*) 'lalim ',lalim(igout) |
---|
| 423 | write(lunout1,*) ' ig l alim_star entr_star detr_star f_star ' |
---|
| 424 | write(lunout1,'(i6,i4,e15.5)') (igout,l,alim_star(igout,l) & |
---|
| 425 | & ,l=1,lalim(igout)+4) |
---|
| 426 | endif |
---|
| 427 | |
---|
| 428 | |
---|
| 429 | |
---|
| 430 | |
---|
| 431 | ! Choix de la fonction d'alimentation utilisee pour la fermeture. |
---|
| 432 | ! Apparemment sans importance |
---|
| 433 | alim_star_clos(:,:)=alim_star(:,:) |
---|
| 434 | alim_star_clos(:,:)=entr_star(:,:)+alim_star(:,:) |
---|
[1998] | 435 | ! |
---|
| 436 | !CR Appel de la fermeture seche |
---|
| 437 | if (iflag_thermals_closure.eq.1) then |
---|
[1403] | 438 | |
---|
[4094] | 439 | CALL thermcell_closure(ngrid,nlay,r_aspect_thermals,ptimestep,rho, & |
---|
| 440 | & zlev,lalim,alim_star_clos,zmax_sec,wmax_sec,f) |
---|
[878] | 441 | |
---|
[1403] | 442 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 443 | ! Appel avec les zmax et wmax tenant compte de la condensation |
---|
| 444 | ! Semble moins bien marcher |
---|
[1998] | 445 | else if (iflag_thermals_closure.eq.2) then |
---|
| 446 | |
---|
| 447 | CALL thermcell_closure(ngrid,nlay,r_aspect_thermals,ptimestep,rho, & |
---|
[4094] | 448 | & zlev,lalim,alim_star,zmax,wmax,f) |
---|
[1998] | 449 | |
---|
[4094] | 450 | |
---|
[1998] | 451 | endif |
---|
| 452 | |
---|
[1403] | 453 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 454 | |
---|
[938] | 455 | if(prt_level.ge.1)print*,'thermcell_closure apres thermcell_closure' |
---|
[878] | 456 | |
---|
[972] | 457 | if (tau_thermals>1.) then |
---|
| 458 | lambda=exp(-ptimestep/tau_thermals) |
---|
| 459 | f0=(1.-lambda)*f+lambda*f0 |
---|
| 460 | else |
---|
| 461 | f0=f |
---|
| 462 | endif |
---|
| 463 | |
---|
| 464 | ! Test valable seulement en 1D mais pas genant |
---|
| 465 | if (.not. (f0(1).ge.0.) ) then |
---|
[1403] | 466 | abort_message = '.not. (f0(1).ge.0.)' |
---|
[2311] | 467 | CALL abort_physic (modname,abort_message,1) |
---|
[972] | 468 | endif |
---|
| 469 | |
---|
[878] | 470 | !------------------------------------------------------------------------------- |
---|
| 471 | !deduction des flux |
---|
| 472 | |
---|
[972] | 473 | CALL thermcell_flux2(ngrid,nlay,ptimestep,masse, & |
---|
[878] | 474 | & lalim,lmax,alim_star, & |
---|
| 475 | & entr_star,detr_star,f,rhobarz,zlev,zw2,fm,entr, & |
---|
[972] | 476 | & detr,zqla,lev_out,lunout1,igout) |
---|
[4318] | 477 | |
---|
[972] | 478 | !IM 060508 & detr,zqla,zmax,lev_out,lunout,igout) |
---|
[878] | 479 | |
---|
[938] | 480 | if (prt_level.ge.1) print*,'thermcell_main apres thermcell_flux' |
---|
[4678] | 481 | call test_ltherm(ngrid,nlay,pplay,lalim,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_flux lalim ') |
---|
| 482 | call test_ltherm(ngrid,nlay,pplay,lmax ,ztv,p_o,ztva,zqla,f_star,zw2,'thermcell_flux lmax ') |
---|
[878] | 483 | |
---|
| 484 | !------------------------------------------------------------------ |
---|
[972] | 485 | ! On ne prend pas directement les profils issus des calculs precedents |
---|
| 486 | ! mais on s'autorise genereusement une relaxation vers ceci avec |
---|
| 487 | ! une constante de temps tau_thermals (typiquement 1800s). |
---|
| 488 | !------------------------------------------------------------------ |
---|
[878] | 489 | |
---|
[972] | 490 | if (tau_thermals>1.) then |
---|
| 491 | lambda=exp(-ptimestep/tau_thermals) |
---|
| 492 | fm0=(1.-lambda)*fm+lambda*fm0 |
---|
| 493 | entr0=(1.-lambda)*entr+lambda*entr0 |
---|
[1403] | 494 | detr0=(1.-lambda)*detr+lambda*detr0 |
---|
[878] | 495 | else |
---|
| 496 | fm0=fm |
---|
| 497 | entr0=entr |
---|
| 498 | detr0=detr |
---|
| 499 | endif |
---|
| 500 | |
---|
[4396] | 501 | !------------------------------------------------------------------ |
---|
| 502 | ! Calcul de la fraction de l'ascendance |
---|
| 503 | !------------------------------------------------------------------ |
---|
| 504 | do ig=1,ngrid |
---|
| 505 | fraca(ig,1)=0. |
---|
| 506 | fraca(ig,nlay+1)=0. |
---|
| 507 | enddo |
---|
| 508 | do l=2,nlay |
---|
| 509 | do ig=1,ngrid |
---|
| 510 | if (zw2(ig,l).gt.1.e-10) then |
---|
| 511 | fraca(ig,l)=fm(ig,l)/(rhobarz(ig,l)*zw2(ig,l)) |
---|
| 512 | else |
---|
| 513 | fraca(ig,l)=0. |
---|
| 514 | endif |
---|
| 515 | enddo |
---|
| 516 | enddo |
---|
| 517 | |
---|
[972] | 518 | !c------------------------------------------------------------------ |
---|
| 519 | ! calcul du transport vertical |
---|
| 520 | !------------------------------------------------------------------ |
---|
[4381] | 521 | IF (iflag_thermals_down .GT. 0) THEN |
---|
[4438] | 522 | if (debut) print*,'WARNING !!! routine thermcell_down en cours de developpement' |
---|
[4441] | 523 | entrdn=fact_thermals_down*detr0 |
---|
| 524 | detrdn=fact_thermals_down*entr0 |
---|
[4413] | 525 | ! we want to transport potential temperature, total water and momentum |
---|
| 526 | CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,zthl,zdthladj) |
---|
[4678] | 527 | CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,p_o,pdoadj) |
---|
[4413] | 528 | CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,zu,pduadj) |
---|
| 529 | CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,zv,pdvadj) |
---|
[4396] | 530 | ELSE |
---|
[4381] | 531 | !-------------------------------------------------------------- |
---|
[972] | 532 | |
---|
[4396] | 533 | call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse, & |
---|
| 534 | & zthl,zdthladj,zta,lev_out) |
---|
[4678] | 535 | |
---|
| 536 | do ll=1,nlay |
---|
| 537 | do ig=1,ngrid |
---|
| 538 | z_o(ig,ll)=p_o(ig,ll) |
---|
| 539 | enddo |
---|
| 540 | enddo |
---|
[4396] | 541 | call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse, & |
---|
[4678] | 542 | & z_o,pdoadj,z_oa,lev_out) |
---|
[878] | 543 | |
---|
[4089] | 544 | #ifdef ISO |
---|
[4143] | 545 | ! C Risi: on utilise directement la meme routine |
---|
| 546 | do ixt=1,ntiso |
---|
[4089] | 547 | do ll=1,nlay |
---|
| 548 | DO ig=1,ngrid |
---|
| 549 | xtpo_tmp(ig,ll)=xtpo(ixt,ig,ll) |
---|
| 550 | xtzo_tmp(ig,ll)=xtzo(ixt,ig,ll) |
---|
| 551 | enddo |
---|
| 552 | enddo |
---|
| 553 | call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse, & |
---|
| 554 | & xtpo_tmp,xtpdoadj_tmp,xtzo_tmp,lev_out) |
---|
| 555 | do ll=1,nlay |
---|
| 556 | DO ig=1,ngrid |
---|
| 557 | xtpdoadj(ixt,ig,ll)=xtpdoadj_tmp(ig,ll) |
---|
| 558 | enddo |
---|
| 559 | enddo |
---|
[4143] | 560 | enddo |
---|
[4089] | 561 | #endif |
---|
| 562 | |
---|
| 563 | #ifdef ISO |
---|
| 564 | #ifdef ISOVERIF |
---|
| 565 | DO ll=1,nlay |
---|
| 566 | DO ig=1,ngrid |
---|
| 567 | if (iso_eau.gt.0) then |
---|
| 568 | call iso_verif_egalite(xtpo(iso_eau,ig,ll), & |
---|
[4680] | 569 | & p_o(ig,ll),'thermcell_main 594') |
---|
[4089] | 570 | call iso_verif_egalite(xtpdoadj(iso_eau,ig,ll), & |
---|
| 571 | & pdoadj(ig,ll),'thermcell_main 596') |
---|
| 572 | endif |
---|
| 573 | if (iso_HDO.gt.0) then |
---|
| 574 | call iso_verif_aberrant_encadre(xtpo(iso_hdo,ig,ll) & |
---|
[4678] | 575 | & /p_o(ig,ll),'thermcell_main 610') |
---|
[4089] | 576 | endif |
---|
| 577 | enddo |
---|
| 578 | enddo !DO ll=1,nlay |
---|
| 579 | write(*,*) 'thermcell_main 600 tmp: apres thermcell_dq' |
---|
| 580 | #endif |
---|
| 581 | #endif |
---|
| 582 | |
---|
| 583 | |
---|
[883] | 584 | !------------------------------------------------------------------ |
---|
| 585 | ! calcul du transport vertical du moment horizontal |
---|
| 586 | !------------------------------------------------------------------ |
---|
[878] | 587 | |
---|
[972] | 588 | !IM 090508 |
---|
[1738] | 589 | if (dvdq == 0 ) then |
---|
[883] | 590 | |
---|
[878] | 591 | ! Calcul du transport de V tenant compte d'echange par gradient |
---|
| 592 | ! de pression horizontal avec l'environnement |
---|
| 593 | |
---|
| 594 | call thermcell_dv2(ngrid,nlay,ptimestep,fm0,entr0,masse & |
---|
[1738] | 595 | ! & ,fraca*dvdq,zmax & |
---|
| 596 | & ,fraca,zmax & |
---|
[972] | 597 | & ,zu,zv,pduadj,pdvadj,zua,zva,lev_out) |
---|
[1403] | 598 | |
---|
[878] | 599 | else |
---|
| 600 | |
---|
| 601 | ! calcul purement conservatif pour le transport de V |
---|
[1738] | 602 | call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse & |
---|
[878] | 603 | & ,zu,pduadj,zua,lev_out) |
---|
[1738] | 604 | call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse & |
---|
[878] | 605 | & ,zv,pdvadj,zva,lev_out) |
---|
[1738] | 606 | |
---|
[878] | 607 | endif |
---|
[4396] | 608 | ENDIF |
---|
[878] | 609 | |
---|
| 610 | ! print*,'13 OK convect8' |
---|
| 611 | do l=1,nlay |
---|
| 612 | do ig=1,ngrid |
---|
| 613 | pdtadj(ig,l)=zdthladj(ig,l)*zpspsk(ig,l) |
---|
| 614 | enddo |
---|
| 615 | enddo |
---|
| 616 | |
---|
[972] | 617 | if (prt_level.ge.1) print*,'14 OK convect8' |
---|
[878] | 618 | !------------------------------------------------------------------ |
---|
| 619 | ! Calculs de diagnostiques pour les sorties |
---|
| 620 | !------------------------------------------------------------------ |
---|
| 621 | !calcul de fraca pour les sorties |
---|
| 622 | |
---|
| 623 | if (sorties) then |
---|
[972] | 624 | if (prt_level.ge.1) print*,'14a OK convect8' |
---|
[878] | 625 | ! calcul du niveau de condensation |
---|
| 626 | ! initialisation |
---|
| 627 | do ig=1,ngrid |
---|
[879] | 628 | nivcon(ig)=0 |
---|
[878] | 629 | zcon(ig)=0. |
---|
| 630 | enddo |
---|
| 631 | !nouveau calcul |
---|
| 632 | do ig=1,ngrid |
---|
[4678] | 633 | CHI=zh(ig,1)/(1669.0-122.0*z_o(ig,1)/zqsat(ig,1)-zh(ig,1)) |
---|
| 634 | pcon(ig)=pplay(ig,1)*(z_o(ig,1)/zqsat(ig,1))**CHI |
---|
[878] | 635 | enddo |
---|
[1403] | 636 | !IM do k=1,nlay |
---|
| 637 | do k=1,nlay-1 |
---|
[878] | 638 | do ig=1,ngrid |
---|
| 639 | if ((pcon(ig).le.pplay(ig,k)) & |
---|
| 640 | & .and.(pcon(ig).gt.pplay(ig,k+1))) then |
---|
| 641 | zcon2(ig)=zlay(ig,k)-(pcon(ig)-pplay(ig,k))/(RG*rho(ig,k))/100. |
---|
| 642 | endif |
---|
| 643 | enddo |
---|
| 644 | enddo |
---|
[1403] | 645 | !IM |
---|
[1494] | 646 | ierr=0 |
---|
[1403] | 647 | do ig=1,ngrid |
---|
| 648 | if (pcon(ig).le.pplay(ig,nlay)) then |
---|
| 649 | zcon2(ig)=zlay(ig,nlay)-(pcon(ig)-pplay(ig,nlay))/(RG*rho(ig,nlay))/100. |
---|
[1494] | 650 | ierr=1 |
---|
| 651 | endif |
---|
| 652 | enddo |
---|
| 653 | if (ierr==1) then |
---|
[1403] | 654 | abort_message = 'thermcellV0_main: les thermiques vont trop haut ' |
---|
[2311] | 655 | CALL abort_physic (modname,abort_message,1) |
---|
[1494] | 656 | endif |
---|
| 657 | |
---|
[972] | 658 | if (prt_level.ge.1) print*,'14b OK convect8' |
---|
[878] | 659 | do k=nlay,1,-1 |
---|
| 660 | do ig=1,ngrid |
---|
| 661 | if (zqla(ig,k).gt.1e-10) then |
---|
| 662 | nivcon(ig)=k |
---|
| 663 | zcon(ig)=zlev(ig,k) |
---|
| 664 | endif |
---|
| 665 | enddo |
---|
| 666 | enddo |
---|
[972] | 667 | if (prt_level.ge.1) print*,'14c OK convect8' |
---|
[878] | 668 | !calcul des moments |
---|
| 669 | !initialisation |
---|
| 670 | do l=1,nlay |
---|
| 671 | do ig=1,ngrid |
---|
| 672 | q2(ig,l)=0. |
---|
| 673 | wth2(ig,l)=0. |
---|
| 674 | wth3(ig,l)=0. |
---|
| 675 | ratqscth(ig,l)=0. |
---|
| 676 | ratqsdiff(ig,l)=0. |
---|
| 677 | enddo |
---|
| 678 | enddo |
---|
[972] | 679 | if (prt_level.ge.1) print*,'14d OK convect8' |
---|
[1146] | 680 | if (prt_level.ge.10)write(lunout,*) & |
---|
| 681 | & 'WARNING thermcell_main wth2=0. si zw2 > 1.e-10' |
---|
[878] | 682 | do l=1,nlay |
---|
| 683 | do ig=1,ngrid |
---|
| 684 | zf=fraca(ig,l) |
---|
| 685 | zf2=zf/(1.-zf) |
---|
[972] | 686 | ! |
---|
[1403] | 687 | thetath2(ig,l)=zf2*(ztla(ig,l)-zthl(ig,l))**2 |
---|
[972] | 688 | if(zw2(ig,l).gt.1.e-10) then |
---|
| 689 | wth2(ig,l)=zf2*(zw2(ig,l))**2 |
---|
| 690 | else |
---|
| 691 | wth2(ig,l)=0. |
---|
| 692 | endif |
---|
[878] | 693 | wth3(ig,l)=zf2*(1-2.*fraca(ig,l))/(1-fraca(ig,l)) & |
---|
| 694 | & *zw2(ig,l)*zw2(ig,l)*zw2(ig,l) |
---|
[4678] | 695 | q2(ig,l)=zf2*(zqta(ig,l)*1000.-p_o(ig,l)*1000.)**2 |
---|
| 696 | !test: on calcul q2/p_o=ratqsc |
---|
| 697 | ratqscth(ig,l)=sqrt(max(q2(ig,l),1.e-6)/(p_o(ig,l)*1000.)) |
---|
[878] | 698 | enddo |
---|
| 699 | enddo |
---|
[1403] | 700 | !calcul des flux: q, thetal et thetav |
---|
| 701 | do l=1,nlay |
---|
| 702 | do ig=1,ngrid |
---|
[4678] | 703 | wq(ig,l)=fraca(ig,l)*zw2(ig,l)*(zqta(ig,l)*1000.-p_o(ig,l)*1000.) |
---|
[1403] | 704 | wthl(ig,l)=fraca(ig,l)*zw2(ig,l)*(ztla(ig,l)-zthl(ig,l)) |
---|
| 705 | wthv(ig,l)=fraca(ig,l)*zw2(ig,l)*(ztva(ig,l)-ztv(ig,l)) |
---|
| 706 | enddo |
---|
[879] | 707 | enddo |
---|
[1638] | 708 | |
---|
[878] | 709 | !calcul du ratqscdiff |
---|
[972] | 710 | if (prt_level.ge.1) print*,'14e OK convect8' |
---|
[878] | 711 | var=0. |
---|
| 712 | vardiff=0. |
---|
| 713 | ratqsdiff(:,:)=0. |
---|
[1494] | 714 | |
---|
[4089] | 715 | do l=1,nlay |
---|
[1494] | 716 | do ig=1,ngrid |
---|
| 717 | if (l<=lalim(ig)) then |
---|
[878] | 718 | var=var+alim_star(ig,l)*zqta(ig,l)*1000. |
---|
[1494] | 719 | endif |
---|
[878] | 720 | enddo |
---|
| 721 | enddo |
---|
[1494] | 722 | |
---|
[972] | 723 | if (prt_level.ge.1) print*,'14f OK convect8' |
---|
[1494] | 724 | |
---|
[4089] | 725 | do l=1,nlay |
---|
[1494] | 726 | do ig=1,ngrid |
---|
| 727 | if (l<=lalim(ig)) then |
---|
| 728 | zf=fraca(ig,l) |
---|
| 729 | zf2=zf/(1.-zf) |
---|
| 730 | vardiff=vardiff+alim_star(ig,l)*(zqta(ig,l)*1000.-var)**2 |
---|
| 731 | endif |
---|
| 732 | enddo |
---|
[878] | 733 | enddo |
---|
[1494] | 734 | |
---|
[972] | 735 | if (prt_level.ge.1) print*,'14g OK convect8' |
---|
[4089] | 736 | do l=1,nlay |
---|
| 737 | do ig=1,ngrid |
---|
[4678] | 738 | ratqsdiff(ig,l)=sqrt(vardiff)/(p_o(ig,l)*1000.) |
---|
[4089] | 739 | enddo |
---|
| 740 | enddo |
---|
[878] | 741 | endif |
---|
| 742 | |
---|
[938] | 743 | if (prt_level.ge.1) print*,'thermcell_main FIN OK' |
---|
[878] | 744 | |
---|
[4094] | 745 | RETURN |
---|
[4089] | 746 | end subroutine thermcell_main |
---|
[878] | 747 | |
---|
[4089] | 748 | !============================================================================= |
---|
| 749 | !///////////////////////////////////////////////////////////////////////////// |
---|
| 750 | !============================================================================= |
---|
[4678] | 751 | subroutine test_ltherm(ngrid,nlay,pplay,long,ztv,p_o,ztva, & ! in |
---|
[4089] | 752 | & zqla,f_star,zw2,comment) ! in |
---|
| 753 | !============================================================================= |
---|
[4590] | 754 | USE lmdz_thermcell_ini, ONLY: prt_level |
---|
[938] | 755 | IMPLICIT NONE |
---|
[878] | 756 | |
---|
[4089] | 757 | integer i, k, ngrid,nlay |
---|
[4678] | 758 | real, intent(in), dimension(ngrid,nlay) :: pplay,ztv,p_o,ztva,zqla |
---|
[4089] | 759 | real, intent(in), dimension(ngrid,nlay) :: f_star,zw2 |
---|
| 760 | integer, intent(in), dimension(ngrid) :: long |
---|
[878] | 761 | real seuil |
---|
| 762 | character*21 comment |
---|
[4133] | 763 | |
---|
[4089] | 764 | seuil=0.25 |
---|
[878] | 765 | |
---|
[938] | 766 | if (prt_level.ge.1) THEN |
---|
| 767 | print*,'WARNING !!! TEST ',comment |
---|
| 768 | endif |
---|
[879] | 769 | return |
---|
| 770 | |
---|
[878] | 771 | ! test sur la hauteur des thermiques ... |
---|
[4089] | 772 | do i=1,ngrid |
---|
[972] | 773 | !IMtemp if (pplay(i,long(i)).lt.seuil*pplev(i,1)) then |
---|
| 774 | if (prt_level.ge.10) then |
---|
[878] | 775 | print*,'WARNING ',comment,' au point ',i,' K= ',long(i) |
---|
| 776 | print*,' K P(MB) THV(K) Qenv(g/kg)THVA QLA(g/kg) F* W2' |
---|
[4089] | 777 | do k=1,nlay |
---|
[4678] | 778 | write(6,'(i3,7f10.3)') k,pplay(i,k),ztv(i,k),1000*p_o(i,k),ztva(i,k),1000*zqla(i,k),f_star(i,k),zw2(i,k) |
---|
[878] | 779 | enddo |
---|
[972] | 780 | endif |
---|
[878] | 781 | enddo |
---|
| 782 | |
---|
| 783 | |
---|
| 784 | return |
---|
| 785 | end |
---|
| 786 | |
---|
[4089] | 787 | ! nrlmd le 10/04/2012 Transport de la TKE par le thermique moyen pour la fermeture en ALP |
---|
| 788 | ! On transporte pbl_tke pour donner therm_tke |
---|
[4143] | 789 | ! Copie conforme de la subroutine DTKE dans physiq.F ecrite par Frederic Hourdin |
---|
[4089] | 790 | |
---|
| 791 | !======================================================================= |
---|
| 792 | !/////////////////////////////////////////////////////////////////////// |
---|
| 793 | !======================================================================= |
---|
| 794 | |
---|
| 795 | subroutine thermcell_tke_transport( & |
---|
| 796 | & ngrid,nlay,ptimestep,fm0,entr0,rg,pplev, & ! in |
---|
| 797 | & therm_tke_max) ! out |
---|
[4590] | 798 | USE lmdz_thermcell_ini, ONLY: prt_level |
---|
[1638] | 799 | implicit none |
---|
| 800 | |
---|
| 801 | !======================================================================= |
---|
| 802 | ! |
---|
| 803 | ! Calcul du transport verticale dans la couche limite en presence |
---|
| 804 | ! de "thermiques" explicitement representes |
---|
| 805 | ! calcul du dq/dt une fois qu'on connait les ascendances |
---|
| 806 | ! |
---|
| 807 | !======================================================================= |
---|
| 808 | |
---|
[4089] | 809 | integer ngrid,nlay |
---|
[1638] | 810 | |
---|
[4089] | 811 | real, intent(in) :: ptimestep |
---|
| 812 | real, intent(in), dimension(ngrid,nlay+1) :: fm0,pplev |
---|
| 813 | real, intent(in), dimension(ngrid,nlay) :: entr0 |
---|
| 814 | real, intent(in) :: rg |
---|
| 815 | real, intent(out), dimension(ngrid,nlay) :: therm_tke_max |
---|
| 816 | |
---|
[1638] | 817 | real detr0(ngrid,nlay) |
---|
[4089] | 818 | real masse0(ngrid,nlay) |
---|
[1638] | 819 | real masse(ngrid,nlay),fm(ngrid,nlay+1) |
---|
| 820 | real entr(ngrid,nlay) |
---|
| 821 | real q(ngrid,nlay) |
---|
| 822 | integer lev_out ! niveau pour les print |
---|
| 823 | |
---|
| 824 | real qa(ngrid,nlay),detr(ngrid,nlay),wqd(ngrid,nlay+1) |
---|
| 825 | integer ig,k |
---|
| 826 | |
---|
| 827 | |
---|
| 828 | lev_out=0 |
---|
| 829 | |
---|
| 830 | |
---|
| 831 | if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0' |
---|
| 832 | |
---|
| 833 | ! calcul du detrainement |
---|
| 834 | do k=1,nlay |
---|
| 835 | detr0(:,k)=fm0(:,k)-fm0(:,k+1)+entr0(:,k) |
---|
| 836 | masse0(:,k)=(pplev(:,k)-pplev(:,k+1))/RG |
---|
| 837 | enddo |
---|
| 838 | |
---|
| 839 | |
---|
| 840 | ! Decalage vertical des entrainements et detrainements. |
---|
| 841 | masse(:,1)=0.5*masse0(:,1) |
---|
| 842 | entr(:,1)=0.5*entr0(:,1) |
---|
| 843 | detr(:,1)=0.5*detr0(:,1) |
---|
| 844 | fm(:,1)=0. |
---|
| 845 | do k=1,nlay-1 |
---|
| 846 | masse(:,k+1)=0.5*(masse0(:,k)+masse0(:,k+1)) |
---|
| 847 | entr(:,k+1)=0.5*(entr0(:,k)+entr0(:,k+1)) |
---|
| 848 | detr(:,k+1)=0.5*(detr0(:,k)+detr0(:,k+1)) |
---|
| 849 | fm(:,k+1)=fm(:,k)+entr(:,k)-detr(:,k) |
---|
| 850 | enddo |
---|
| 851 | fm(:,nlay+1)=0. |
---|
| 852 | |
---|
| 853 | |
---|
| 854 | q(:,:)=therm_tke_max(:,:) |
---|
| 855 | !!! nrlmd le 16/09/2010 |
---|
| 856 | do ig=1,ngrid |
---|
| 857 | qa(ig,1)=q(ig,1) |
---|
| 858 | enddo |
---|
| 859 | !!! |
---|
| 860 | |
---|
| 861 | if (1==1) then |
---|
| 862 | do k=2,nlay |
---|
| 863 | do ig=1,ngrid |
---|
| 864 | if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt. & |
---|
| 865 | & 1.e-5*masse(ig,k)) then |
---|
| 866 | qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k)) & |
---|
| 867 | & /(fm(ig,k+1)+detr(ig,k)) |
---|
| 868 | else |
---|
| 869 | qa(ig,k)=q(ig,k) |
---|
| 870 | endif |
---|
| 871 | if (qa(ig,k).lt.0.) then |
---|
| 872 | ! print*,'qa<0!!!' |
---|
| 873 | endif |
---|
| 874 | if (q(ig,k).lt.0.) then |
---|
| 875 | ! print*,'q<0!!!' |
---|
| 876 | endif |
---|
| 877 | enddo |
---|
| 878 | enddo |
---|
| 879 | |
---|
| 880 | ! Calcul du flux subsident |
---|
| 881 | |
---|
| 882 | do k=2,nlay |
---|
| 883 | do ig=1,ngrid |
---|
| 884 | wqd(ig,k)=fm(ig,k)*q(ig,k) |
---|
| 885 | if (wqd(ig,k).lt.0.) then |
---|
| 886 | ! print*,'wqd<0!!!' |
---|
| 887 | endif |
---|
| 888 | enddo |
---|
| 889 | enddo |
---|
| 890 | do ig=1,ngrid |
---|
| 891 | wqd(ig,1)=0. |
---|
| 892 | wqd(ig,nlay+1)=0. |
---|
| 893 | enddo |
---|
| 894 | |
---|
| 895 | ! Calcul des tendances |
---|
| 896 | do k=1,nlay |
---|
| 897 | do ig=1,ngrid |
---|
| 898 | q(ig,k)=q(ig,k)+(detr(ig,k)*qa(ig,k)-entr(ig,k)*q(ig,k) & |
---|
| 899 | & -wqd(ig,k)+wqd(ig,k+1)) & |
---|
| 900 | & *ptimestep/masse(ig,k) |
---|
| 901 | enddo |
---|
| 902 | enddo |
---|
| 903 | |
---|
| 904 | endif |
---|
| 905 | |
---|
| 906 | therm_tke_max(:,:)=q(:,:) |
---|
| 907 | |
---|
| 908 | return |
---|
| 909 | !!! fin nrlmd le 10/04/2012 |
---|
| 910 | end |
---|
| 911 | |
---|
[4590] | 912 | END MODULE lmdz_thermcell_main |
---|