Directory: | ./ |
---|---|
File: | phys/thermcell_alp.f90 |
Date: | 2022-01-11 19:19:34 |
Exec | Total | Coverage | |
---|---|---|---|
Lines: | 144 | 152 | 94.7% |
Branches: | 123 | 142 | 86.6% |
Line | Branch | Exec | Source |
---|---|---|---|
1 | ! $Id: thermcell_main.F90 2351 2015-08-25 15:14:59Z emillour $ | ||
2 | ! | ||
3 | 18973279 | SUBROUTINE thermcell_alp(ngrid,nlay,ptimestep & | |
4 | 480 | & ,pplay,pplev & | |
5 | & ,fm0,entr0,lmax & | ||
6 | 480 | & ,ale_bl,alp_bl,lalim_conv,wght_th & | |
7 | & ,zw2,fraca & | ||
8 | !!! ncessaire en plus | ||
9 | & ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax & | ||
10 | !!! nrlmd le 10/04/2012 | ||
11 | 480 | & ,pbl_tke,pctsrf,omega,airephy & | |
12 | 480 | & ,zlcl,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 & | |
13 | & ,n2,s2,ale_bl_stat & | ||
14 | 480 | & ,therm_tke_max,env_tke_max & | |
15 | & ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke & | ||
16 | & ,alp_bl_conv,alp_bl_stat & | ||
17 | !!! fin nrlmd le 10/04/2012 | ||
18 | &) | ||
19 | |||
20 | USE dimphy | ||
21 | USE indice_sol_mod | ||
22 | IMPLICIT NONE | ||
23 | |||
24 | !======================================================================= | ||
25 | ! Auteurs: Frederic Hourdin, Catherine Rio, Anne Mathieu | ||
26 | ! Version du 09.02.07 | ||
27 | ! Calcul du transport vertical dans la couche limite en presence | ||
28 | ! de "thermiques" explicitement representes avec processus nuageux | ||
29 | ! | ||
30 | ! Reecriture a partir d'un listing papier a Habas, le 14/02/00 | ||
31 | ! | ||
32 | ! le thermique est suppose homogene et dissipe par melange avec | ||
33 | ! son environnement. la longueur l_mix controle l'efficacite du | ||
34 | ! melange | ||
35 | ! | ||
36 | ! Le calcul du transport des differentes especes se fait en prenant | ||
37 | ! en compte: | ||
38 | ! 1. un flux de masse montant | ||
39 | ! 2. un flux de masse descendant | ||
40 | ! 3. un entrainement | ||
41 | ! 4. un detrainement | ||
42 | ! | ||
43 | ! Modif 2013/01/04 (FH hourdin@lmd.jussieu.fr) | ||
44 | ! Introduction of an implicit computation of vertical advection in | ||
45 | ! the environment of thermal plumes in thermcell_dq | ||
46 | ! impl = 0 : explicit, 1 : implicit, -1 : old version | ||
47 | ! controled by iflag_thermals = | ||
48 | ! 15, 16 run with impl=-1 : numerical convergence with NPv3 | ||
49 | ! 17, 18 run with impl=1 : more stable | ||
50 | ! 15 and 17 correspond to the activation of the stratocumulus "bidouille" | ||
51 | ! | ||
52 | !======================================================================= | ||
53 | !----------------------------------------------------------------------- | ||
54 | ! declarations: | ||
55 | ! ------------- | ||
56 | |||
57 | ! | ||
58 | ! $Header$ | ||
59 | ! | ||
60 | ! ATTENTION!!!!: ce fichier include est compatible format fixe/format libre | ||
61 | ! veillez � n'utiliser que des ! pour les commentaires | ||
62 | ! et � bien positionner les & des lignes de continuation | ||
63 | ! (les placer en colonne 6 et en colonne 73) | ||
64 | ! | ||
65 | ! | ||
66 | ! A1.0 Fundamental constants | ||
67 | REAL RPI,RCLUM,RHPLA,RKBOL,RNAVO | ||
68 | ! A1.1 Astronomical constants | ||
69 | REAL RDAY,REA,REPSM,RSIYEA,RSIDAY,ROMEGA | ||
70 | ! A1.1.bis Constantes concernant l'orbite de la Terre: | ||
71 | REAL R_ecc, R_peri, R_incl | ||
72 | ! A1.2 Geoide | ||
73 | REAL RA,RG,R1SA | ||
74 | ! A1.3 Radiation | ||
75 | ! REAL RSIGMA,RI0 | ||
76 | REAL RSIGMA | ||
77 | ! A1.4 Thermodynamic gas phase | ||
78 | REAL RMO3,RMCO2,RMC,RMCH4,RMN2O,RMCFC11,RMCFC12 | ||
79 | REAL R,RMD,RMV,RD,RV,RCPD,RCPV,RCVD,RCVV | ||
80 | REAL RKAPPA,RETV, eps_w | ||
81 | ! A1.5,6 Thermodynamic liquid,solid phases | ||
82 | REAL RCW,RCS | ||
83 | ! A1.7 Thermodynamic transition of phase | ||
84 | REAL RLVTT,RLSTT,RLMLT,RTT,RATM | ||
85 | ! A1.8 Curve of saturation | ||
86 | REAL RESTT,RALPW,RBETW,RGAMW,RALPS,RBETS,RGAMS | ||
87 | REAL RALPD,RBETD,RGAMD | ||
88 | ! | ||
89 | COMMON/YOMCST/RPI ,RCLUM ,RHPLA ,RKBOL ,RNAVO & | ||
90 | & ,RDAY ,REA ,REPSM ,RSIYEA,RSIDAY,ROMEGA & | ||
91 | & ,R_ecc, R_peri, R_incl & | ||
92 | & ,RA ,RG ,R1SA & | ||
93 | & ,RSIGMA & | ||
94 | & ,R ,RMD ,RMV ,RD ,RV ,RCPD & | ||
95 | & ,RMO3 ,RMCO2 ,RMC ,RMCH4 ,RMN2O ,RMCFC11 ,RMCFC12 & | ||
96 | & ,RCPV ,RCVD ,RCVV ,RKAPPA,RETV, eps_w & | ||
97 | & ,RCW ,RCS & | ||
98 | & ,RLVTT ,RLSTT ,RLMLT ,RTT ,RATM & | ||
99 | & ,RESTT ,RALPW ,RBETW ,RGAMW ,RALPS ,RBETS ,RGAMS & | ||
100 | & ,RALPD ,RBETD ,RGAMD | ||
101 | ! ------------------------------------------------------------------ | ||
102 | !$OMP THREADPRIVATE(/YOMCST/) | ||
103 | ! | ||
104 | ! $Id: YOETHF.h 2799 2017-02-24 18:50:33Z jyg $ | ||
105 | ! | ||
106 | ! ATTENTION!!!!: ce fichier include est compatible format fixe/format libre | ||
107 | ! veillez n'utiliser que des ! pour les commentaires | ||
108 | ! et bien positionner les & des lignes de continuation | ||
109 | ! (les placer en colonne 6 et en colonne 73) | ||
110 | ! | ||
111 | !* COMMON *YOETHF* DERIVED CONSTANTS SPECIFIC TO ECMWF THERMODYNAMICS | ||
112 | ! | ||
113 | ! *R__ES* *CONSTANTS USED FOR COMPUTATION OF SATURATION | ||
114 | ! MIXING RATIO OVER LIQUID WATER(*R_LES*) OR | ||
115 | ! ICE(*R_IES*). | ||
116 | ! *RVTMP2* *RVTMP2=RCPV/RCPD-1. | ||
117 | ! *RHOH2O* *DENSITY OF LIQUID WATER. (RATM/100.) | ||
118 | ! | ||
119 | REAL R2ES, R3LES, R3IES, R4LES, R4IES, R5LES, R5IES | ||
120 | REAL RVTMP2, RHOH2O | ||
121 | REAL R5ALVCP,R5ALSCP,RALVDCP,RALSDCP,RALFDCP,RTWAT,RTBER,RTBERCU | ||
122 | REAL RTICE,RTICECU,RTWAT_RTICE_R,RTWAT_RTICECU_R,RKOOP1,RKOOP2 | ||
123 | LOGICAL OK_BAD_ECMWF_THERMO ! If TRUE, then variables set by rrtm/suphec.F90 | ||
124 | ! If FALSE, then variables set by suphel.F90 | ||
125 | COMMON /YOETHF/R2ES, R3LES, R3IES, R4LES, R4IES, R5LES, R5IES, & | ||
126 | & RVTMP2, RHOH2O, & | ||
127 | & R5ALVCP,R5ALSCP,RALVDCP,RALSDCP, & | ||
128 | & RALFDCP,RTWAT,RTBER,RTBERCU, & | ||
129 | & RTICE,RTICECU,RTWAT_RTICE_R,RTWAT_RTICECU_R,RKOOP1,& | ||
130 | & RKOOP2, & | ||
131 | & OK_BAD_ECMWF_THERMO | ||
132 | |||
133 | !$OMP THREADPRIVATE(/YOETHF/) | ||
134 | ! | ||
135 | ! $Header$ | ||
136 | ! | ||
137 | ! | ||
138 | ! ATTENTION!!!!: ce fichier include est compatible format fixe/format libre | ||
139 | ! veillez n'utiliser que des ! pour les commentaires | ||
140 | ! et bien positionner les & des lignes de continuation | ||
141 | ! (les placer en colonne 6 et en colonne 73) | ||
142 | ! | ||
143 | ! ------------------------------------------------------------------ | ||
144 | ! This COMDECK includes the Thermodynamical functions for the cy39 | ||
145 | ! ECMWF Physics package. | ||
146 | ! Consistent with YOMCST Basic physics constants, assuming the | ||
147 | ! partial pressure of water vapour is given by a first order | ||
148 | ! Taylor expansion of Qs(T) w.r.t. to Temperature, using constants | ||
149 | ! in YOETHF | ||
150 | ! ------------------------------------------------------------------ | ||
151 | REAL PTARG, PDELARG, P5ARG, PQSARG, PCOARG | ||
152 | REAL FOEEW, FOEDE, qsats, qsatl, dqsats, dqsatl | ||
153 | LOGICAL thermcep | ||
154 | PARAMETER (thermcep=.TRUE.) | ||
155 | ! | ||
156 | FOEEW ( PTARG,PDELARG ) = EXP ( & | ||
157 | & (R3LES*(1.-PDELARG)+R3IES*PDELARG) * (PTARG-RTT) & | ||
158 | & / (PTARG-(R4LES*(1.-PDELARG)+R4IES*PDELARG)) ) | ||
159 | ! | ||
160 | FOEDE ( PTARG,PDELARG,P5ARG,PQSARG,PCOARG ) = PQSARG*PCOARG*P5ARG & | ||
161 | & / (PTARG-(R4LES*(1.-PDELARG)+R4IES*PDELARG))**2 | ||
162 | ! | ||
163 | qsats(ptarg) = 100.0 * 0.622 * 10.0 & | ||
164 | & ** (2.07023 - 0.00320991 * ptarg & | ||
165 | & - 2484.896 / ptarg + 3.56654 * LOG10(ptarg)) | ||
166 | qsatl(ptarg) = 100.0 * 0.622 * 10.0 & | ||
167 | & ** (23.8319 - 2948.964 / ptarg & | ||
168 | & - 5.028 * LOG10(ptarg) & | ||
169 | & - 29810.16 * EXP( - 0.0699382 * ptarg) & | ||
170 | & + 25.21935 * EXP( - 2999.924 / ptarg)) | ||
171 | ! | ||
172 | dqsats(ptarg,pqsarg) = RLVTT/RCPD*pqsarg * (3.56654/ptarg & | ||
173 | & +2484.896*LOG(10.)/ptarg**2 & | ||
174 | & -0.00320991*LOG(10.)) | ||
175 | dqsatl(ptarg,pqsarg) = RLVTT/RCPD*pqsarg*LOG(10.)* & | ||
176 | & (2948.964/ptarg**2-5.028/LOG(10.)/ptarg & | ||
177 | & +25.21935*2999.924/ptarg**2*EXP(-2999.924/ptarg) & | ||
178 | & +29810.16*0.0699382*EXP(-0.0699382*ptarg)) | ||
179 | integer :: iflag_thermals,nsplit_thermals | ||
180 | |||
181 | !!! nrlmd le 10/04/2012 | ||
182 | integer :: iflag_trig_bl,iflag_clos_bl | ||
183 | integer :: tau_trig_shallow,tau_trig_deep | ||
184 | real :: s_trig | ||
185 | !!! fin nrlmd le 10/04/2012 | ||
186 | |||
187 | real,parameter :: r_aspect_thermals=2.,l_mix_thermals=30. | ||
188 | real :: alp_bl_k | ||
189 | real :: tau_thermals,fact_thermals_ed_dz | ||
190 | integer,parameter :: w2di_thermals=0 | ||
191 | integer :: isplit | ||
192 | |||
193 | integer :: iflag_coupl,iflag_clos,iflag_wake | ||
194 | integer :: iflag_thermals_ed,iflag_thermals_optflux,iflag_thermals_closure | ||
195 | |||
196 | common/ctherm1/iflag_thermals,nsplit_thermals,iflag_thermals_closure | ||
197 | common/ctherm2/tau_thermals,alp_bl_k,fact_thermals_ed_dz | ||
198 | common/ctherm4/iflag_coupl,iflag_clos,iflag_wake | ||
199 | common/ctherm5/iflag_thermals_ed,iflag_thermals_optflux | ||
200 | |||
201 | !!! nrlmd le 10/04/2012 | ||
202 | common/ctherm6/iflag_trig_bl,iflag_clos_bl | ||
203 | common/ctherm7/tau_trig_shallow,tau_trig_deep | ||
204 | common/ctherm8/s_trig | ||
205 | !!! fin nrlmd le 10/04/2012 | ||
206 | |||
207 | !$OMP THREADPRIVATE(/ctherm1/,/ctherm2/,/ctherm4/,/ctherm5/) | ||
208 | !$OMP THREADPRIVATE(/ctherm6/,/ctherm7/,/ctherm8/) | ||
209 | |||
210 | ! arguments: | ||
211 | ! ---------- | ||
212 | |||
213 | !IM 140508 | ||
214 | |||
215 | INTEGER ngrid,nlay | ||
216 | real ptimestep | ||
217 | REAL pplay(ngrid,nlay),pplev(ngrid,nlay+1) | ||
218 | |||
219 | ! local: | ||
220 | ! ------ | ||
221 | |||
222 | |||
223 | REAL susqr2pi, reuler | ||
224 | |||
225 | INTEGER ig,k,l | ||
226 | INTEGER lmax(klon),lalim(klon) | ||
227 | real zmax(klon),zw2(klon,klev+1) | ||
228 | |||
229 | !on garde le zmax du pas de temps precedent | ||
230 | |||
231 | |||
232 | real fraca(klon,klev+1) | ||
233 | real wth3(klon,klev) | ||
234 | ! FH probleme de dimensionnement avec l'allocation dynamique | ||
235 | ! common/comtherm/thetath2,wth2 | ||
236 | real rhobarz(klon,klev) | ||
237 | |||
238 | real wmax_sec(klon) | ||
239 | real fm0(klon,klev+1),entr0(klon,klev) | ||
240 | real fm(klon,klev+1) | ||
241 | |||
242 | !niveau de condensation | ||
243 | real pcon(klon) | ||
244 | |||
245 | real alim_star(klon,klev) | ||
246 | |||
247 | !!! nrlmd le 10/04/2012 | ||
248 | |||
249 | !------Entr�es | ||
250 | real pbl_tke(klon,klev+1,nbsrf) | ||
251 | real pctsrf(klon,nbsrf) | ||
252 | real omega(klon,klev) | ||
253 | real airephy(klon) | ||
254 | !------Sorties | ||
255 | real zlcl(klon),fraca0(klon),w0(klon),w_conv(klon) | ||
256 | real therm_tke_max0(klon),env_tke_max0(klon) | ||
257 | real n2(klon),s2(klon) | ||
258 | real ale_bl_stat(klon) | ||
259 | real therm_tke_max(klon,klev),env_tke_max(klon,klev) | ||
260 | real alp_bl_det(klon),alp_bl_fluct_m(klon),alp_bl_fluct_tke(klon),alp_bl_conv(klon),alp_bl_stat(klon) | ||
261 | !------Local | ||
262 | integer nsrf | ||
263 | 960 | real rhobarz0(klon) ! Densit� au LCL | |
264 | 960 | logical ok_lcl(klon) ! Existence du LCL des thermiques | |
265 | 960 | integer klcl(klon) ! Niveau du LCL | |
266 | 960 | real interp(klon) ! Coef d'interpolation pour le LCL | |
267 | !--Triggering | ||
268 | real Su ! Surface unit�: celle d'un updraft �l�mentaire | ||
269 | parameter(Su=4e4) | ||
270 | real hcoef ! Coefficient directeur pour le calcul de s2 | ||
271 | parameter(hcoef=1) | ||
272 | real hmincoef ! Coefficient directeur pour l'ordonn�e � l'origine pour le calcul de s2 | ||
273 | parameter(hmincoef=0.3) | ||
274 | real eps1 ! Fraction de surface occup�e par la population 1 : eps1=n1*s1/(fraca0*Sd) | ||
275 | parameter(eps1=0.3) | ||
276 | 960 | real hmin(ngrid) ! Ordonn�e � l'origine pour le calcul de s2 | |
277 | 960 | real zmax_moy(ngrid) ! Hauteur moyenne des thermiques : zmax_moy = zlcl + 0.33 (zmax-zlcl) | |
278 | real zmax_moy_coef | ||
279 | parameter(zmax_moy_coef=0.33) | ||
280 | 960 | real depth(klon) ! Epaisseur moyenne du cumulus | |
281 | 960 | real w_max(klon) ! Vitesse max statistique | |
282 | 960 | real s_max(klon) | |
283 | !--Closure | ||
284 | 960 | real pbl_tke_max(klon,klev) ! Profil de TKE moyenne | |
285 | 960 | real pbl_tke_max0(klon) ! TKE moyenne au LCL | |
286 | 960 | real w_ls(klon,klev) ! Vitesse verticale grande �chelle (m/s) | |
287 | real coef_m ! On consid�re un rendement pour alp_bl_fluct_m | ||
288 | parameter(coef_m=1.) | ||
289 | real coef_tke ! On consid�re un rendement pour alp_bl_fluct_tke | ||
290 | parameter(coef_tke=1.) | ||
291 | |||
292 | !!! fin nrlmd le 10/04/2012 | ||
293 | |||
294 | ! | ||
295 | !nouvelles variables pour la convection | ||
296 | real ale_bl(klon) | ||
297 | real alp_bl(klon) | ||
298 | 960 | real alp_int(klon),dp_int(klon),zdp | |
299 | 480 | real fm_tot(klon) | |
300 | real wght_th(klon,klev) | ||
301 | integer lalim_conv(klon) | ||
302 | !v1d logical therm | ||
303 | !v1d save therm | ||
304 | |||
305 | |||
306 | !------------------------------------------------------------ | ||
307 | ! Initialize output arrays related to stochastic triggering | ||
308 | !------------------------------------------------------------ | ||
309 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | DO ig = 1,klon |
310 | 477120 | zlcl(ig) = 0. | |
311 | 477120 | fraca0(ig) = 0. | |
312 | 477120 | w0(ig) = 0. | |
313 | 477120 | w_conv(ig) = 0. | |
314 | 477120 | therm_tke_max0(ig) = 0. | |
315 | 477120 | env_tke_max0(ig) = 0. | |
316 | 477120 | n2(ig) = 0. | |
317 | 477120 | s2(ig) = 0. | |
318 | 477120 | ale_bl_stat(ig) = 0. | |
319 | 477120 | alp_bl_det(ig) = 0. | |
320 | 477120 | alp_bl_fluct_m(ig) = 0. | |
321 | 477120 | alp_bl_fluct_tke(ig) = 0. | |
322 | 477120 | alp_bl_conv(ig) = 0. | |
323 | 477600 | alp_bl_stat(ig) = 0. | |
324 | ENDDO | ||
325 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | DO l = 1,klev |
326 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | DO ig = 1,klon |
327 | 18607680 | therm_tke_max(ig,l) = 0. | |
328 | 18626400 | env_tke_max(ig,l) = 0. | |
329 | ENDDO | ||
330 | ENDDO | ||
331 | !------------------------------------------------------------ | ||
332 | |||
333 | |||
334 | !------------Test sur le LCL des thermiques | ||
335 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | do ig=1,ngrid |
336 | 477120 | ok_lcl(ig)=.false. | |
337 |
3/4✓ Branch 0 taken 477120 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 362719 times.
✓ Branch 3 taken 114401 times.
|
477600 | if ( (pcon(ig) .gt. pplay(ig,klev-1)) .and. (pcon(ig) .lt. pplay(ig,1)) ) ok_lcl(ig)=.true. |
338 | enddo | ||
339 | |||
340 | !------------Localisation des niveaux entourant le LCL et du coef d'interpolation | ||
341 |
2/2✓ Branch 0 taken 18240 times.
✓ Branch 1 taken 480 times.
|
18720 | do l=1,nlay-1 |
342 |
2/2✓ Branch 0 taken 18130560 times.
✓ Branch 1 taken 18240 times.
|
18149280 | do ig=1,ngrid |
343 |
2/2✓ Branch 0 taken 13783322 times.
✓ Branch 1 taken 4347238 times.
|
18148800 | if (ok_lcl(ig)) then |
344 | !ATTENTION,zw2 calcule en pplev | ||
345 | ! if ((pplay(ig,l) .ge. pcon(ig)) .and. (pplay(ig,l+1) .le. pcon(ig))) then | ||
346 | ! klcl(ig)=l | ||
347 | ! interp(ig)=(pcon(ig)-pplay(ig,klcl(ig)))/(pplay(ig,klcl(ig)+1)-pplay(ig,klcl(ig))) | ||
348 | ! endif | ||
349 |
4/4✓ Branch 0 taken 1599246 times.
✓ Branch 1 taken 12184076 times.
✓ Branch 2 taken 362719 times.
✓ Branch 3 taken 1236527 times.
|
13783322 | if ((pplev(ig,l) .ge. pcon(ig)) .and. (pplev(ig,l+1) .le. pcon(ig))) then |
350 | 362719 | klcl(ig)=l | |
351 | 362719 | interp(ig)=(pcon(ig)-pplev(ig,klcl(ig)))/(pplev(ig,klcl(ig)+1)-pplev(ig,klcl(ig))) | |
352 | endif | ||
353 | endif | ||
354 | enddo | ||
355 | enddo | ||
356 | |||
357 | !------------Hauteur des thermiques | ||
358 | !!jyg le 27/04/2012 | ||
359 | !! do ig =1,ngrid | ||
360 | !! rhobarz0(ig)=rhobarz(ig,klcl(ig))+(rhobarz(ig,klcl(ig)+1) & | ||
361 | !! & -rhobarz(ig,klcl(ig)))*interp(ig) | ||
362 | !! zlcl(ig)=(pplev(ig,1)-pcon(ig))/(rhobarz0(ig)*RG) | ||
363 | !! if ( (.not.ok_lcl(ig)) .or. (zlcl(ig).gt.zmax(ig)) ) zlcl(ig)=zmax(ig) ! Si zclc > zmax alors on pose zlcl = zmax | ||
364 | !! enddo | ||
365 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | do ig =1,ngrid |
366 | !CR:REHABILITATION ZMAX CONTINU | ||
367 |
2/2✓ Branch 0 taken 362719 times.
✓ Branch 1 taken 114401 times.
|
477600 | if (ok_lcl(ig)) then |
368 | rhobarz0(ig)=rhobarz(ig,klcl(ig))+(rhobarz(ig,klcl(ig)+1) & | ||
369 | 362719 | & -rhobarz(ig,klcl(ig)))*interp(ig) | |
370 | 362719 | zlcl(ig)=(pplev(ig,1)-pcon(ig))/(rhobarz0(ig)*RG) | |
371 | 362719 | zlcl(ig)=min(zlcl(ig),zmax(ig)) ! Si zlcl > zmax alors on pose zlcl = zmax | |
372 | else | ||
373 | 114401 | rhobarz0(ig)=0. | |
374 | 114401 | zlcl(ig)=zmax(ig) | |
375 | endif | ||
376 | enddo | ||
377 | !!jyg fin | ||
378 | |||
379 | !------------Calcul des propri�t�s du thermique au LCL | ||
380 |
1/4✓ Branch 0 taken 480 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
480 | IF ( (iflag_trig_bl.ge.1) .or. (iflag_clos_bl.ge.1) ) THEN |
381 | |||
382 | !-----Initialisation de la TKE moyenne | ||
383 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | do l=1,nlay |
384 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | do ig=1,ngrid |
385 | 18626400 | pbl_tke_max(ig,l)=0. | |
386 | enddo | ||
387 | enddo | ||
388 | |||
389 | !-----Calcul de la TKE moyenne | ||
390 |
2/2✓ Branch 0 taken 1920 times.
✓ Branch 1 taken 480 times.
|
2400 | do nsrf=1,nbsrf |
391 |
2/2✓ Branch 0 taken 74880 times.
✓ Branch 1 taken 1920 times.
|
77280 | do l=1,nlay |
392 |
2/2✓ Branch 0 taken 74430720 times.
✓ Branch 1 taken 74880 times.
|
74507520 | do ig=1,ngrid |
393 | 74505600 | pbl_tke_max(ig,l)=pctsrf(ig,nsrf)*pbl_tke(ig,l,nsrf)+pbl_tke_max(ig,l) | |
394 | enddo | ||
395 | enddo | ||
396 | enddo | ||
397 | |||
398 | !-----Initialisations des TKE dans et hors des thermiques | ||
399 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | do l=1,nlay |
400 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | do ig=1,ngrid |
401 | 18607680 | therm_tke_max(ig,l)=pbl_tke_max(ig,l) | |
402 | 18626400 | env_tke_max(ig,l)=pbl_tke_max(ig,l) | |
403 | enddo | ||
404 | enddo | ||
405 | |||
406 | !-----Calcul de la TKE transport�e par les thermiques : therm_tke_max | ||
407 | call thermcell_tke_transport(ngrid,nlay,ptimestep,fm0,entr0, & | ||
408 | 480 | & rg,pplev,therm_tke_max) | |
409 | ! print *,' thermcell_tke_transport -> ' !!jyg | ||
410 | |||
411 | !-----Calcul des profils verticaux de TKE hors thermiques : env_tke_max, et de la vitesse verticale grande �chelle : W_ls | ||
412 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | do l=1,nlay |
413 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | do ig=1,ngrid |
414 | 18607680 | pbl_tke_max(ig,l)=fraca(ig,l)*therm_tke_max(ig,l)+(1.-fraca(ig,l))*env_tke_max(ig,l) ! Recalcul de TKE moyenne apr�s transport de TKE_TH | |
415 | 18607680 | env_tke_max(ig,l)=(pbl_tke_max(ig,l)-fraca(ig,l)*therm_tke_max(ig,l))/(1.-fraca(ig,l)) ! Recalcul de TKE dans l'environnement apr�s transport de TKE_TH | |
416 | 18626400 | w_ls(ig,l)=-1.*omega(ig,l)/(RG*rhobarz(ig,l)) ! Vitesse verticale de grande �chelle | |
417 | enddo | ||
418 | enddo | ||
419 | ! print *,' apres w_ls = ' !!jyg | ||
420 | |||
421 |
2/2✓ Branch 0 taken 480 times.
✓ Branch 1 taken 477120 times.
|
477600 | do ig=1,ngrid |
422 |
2/2✓ Branch 0 taken 362719 times.
✓ Branch 1 taken 114401 times.
|
477600 | if (ok_lcl(ig)) then |
423 | fraca0(ig)=fraca(ig,klcl(ig))+(fraca(ig,klcl(ig)+1) & | ||
424 | 362719 | & -fraca(ig,klcl(ig)))*interp(ig) | |
425 | w0(ig)=zw2(ig,klcl(ig))+(zw2(ig,klcl(ig)+1) & | ||
426 | 362719 | & -zw2(ig,klcl(ig)))*interp(ig) | |
427 | w_conv(ig)=w_ls(ig,klcl(ig))+(w_ls(ig,klcl(ig)+1) & | ||
428 | 362719 | & -w_ls(ig,klcl(ig)))*interp(ig) | |
429 | therm_tke_max0(ig)=therm_tke_max(ig,klcl(ig)) & | ||
430 | 362719 | & +(therm_tke_max(ig,klcl(ig)+1)-therm_tke_max(ig,klcl(ig)))*interp(ig) | |
431 | env_tke_max0(ig)=env_tke_max(ig,klcl(ig))+(env_tke_max(ig,klcl(ig)+1) & | ||
432 | 362719 | & -env_tke_max(ig,klcl(ig)))*interp(ig) | |
433 | pbl_tke_max0(ig)=pbl_tke_max(ig,klcl(ig))+(pbl_tke_max(ig,klcl(ig)+1) & | ||
434 | 362719 | & -pbl_tke_max(ig,klcl(ig)))*interp(ig) | |
435 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 362719 times.
|
362719 | if (therm_tke_max0(ig).ge.20.) therm_tke_max0(ig)=20. |
436 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 362719 times.
|
362719 | if (env_tke_max0(ig).ge.20.) env_tke_max0(ig)=20. |
437 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 362719 times.
|
362719 | if (pbl_tke_max0(ig).ge.20.) pbl_tke_max0(ig)=20. |
438 | else | ||
439 | 114401 | fraca0(ig)=0. | |
440 | 114401 | w0(ig)=0. | |
441 | !!jyg le 27/04/2012 | ||
442 | !! zlcl(ig)=0. | ||
443 | !! | ||
444 | endif | ||
445 | enddo | ||
446 | |||
447 | ENDIF ! IF ( (iflag_trig_bl.ge.1) .or. (iflag_clos_bl.ge.1) ) | ||
448 | ! print *,'ENDIF ( (iflag_trig_bl.ge.1) .or. (iflag_clos_bl.ge.1) ) ' !!jyg | ||
449 | |||
450 | !------------Triggering------------------ | ||
451 |
1/2✓ Branch 0 taken 480 times.
✗ Branch 1 not taken.
|
480 | IF (iflag_trig_bl.ge.1) THEN |
452 | |||
453 | !-----Initialisations | ||
454 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | depth(:)=0. |
455 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | n2(:)=0. |
456 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | s2(:)=100. ! some low value, arbitrary |
457 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | s_max(:)=0. |
458 | |||
459 | !-----Epaisseur du nuage (depth) et d�termination de la queue du spectre de panaches (n2,s2) et du panache le plus gros (s_max) | ||
460 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | do ig=1,ngrid |
461 | 477120 | zmax_moy(ig)=zlcl(ig)+zmax_moy_coef*(zmax(ig)-zlcl(ig)) | |
462 | 477120 | depth(ig)=zmax_moy(ig)-zlcl(ig) | |
463 | 477120 | hmin(ig)=hmincoef*zlcl(ig) | |
464 |
2/2✓ Branch 0 taken 136195 times.
✓ Branch 1 taken 340925 times.
|
477600 | if (depth(ig).ge.10.) then |
465 | 136195 | s2(ig)=(hcoef*depth(ig)+hmin(ig))**2 | |
466 | 136195 | n2(ig)=(1.-eps1)*fraca0(ig)*airephy(ig)/s2(ig) | |
467 | !! | ||
468 | !!jyg le 27/04/2012 | ||
469 | !! s_max(ig)=s2(ig)*log(n2(ig)) | ||
470 | !! if (n2(ig) .lt. 1) s_max(ig)=0. | ||
471 | 136195 | s_max(ig)=s2(ig)*log(max(n2(ig),1.)) | |
472 | !!fin jyg | ||
473 | else | ||
474 | 340925 | n2(ig)=0. | |
475 | 340925 | s_max(ig)=0. | |
476 | endif | ||
477 | enddo | ||
478 | ! print *,'avant Calcul de Wmax ' !!jyg | ||
479 | |||
480 | !-----Calcul de Wmax et ALE_BL_STAT associ�e | ||
481 | !!jyg le 30/04/2012 | ||
482 | !! do ig=1,ngrid | ||
483 | !! if ( (depth(ig).ge.10.) .and. (s_max(ig).gt.1.) ) then | ||
484 | !! w_max(ig)=w0(ig)*(1.+sqrt(2.*log(s_max(ig)/su)-log(2.*3.14)-log(2.*log(s_max(ig)/su)-log(2.*3.14)))) | ||
485 | !! ale_bl_stat(ig)=0.5*w_max(ig)**2 | ||
486 | !! else | ||
487 | !! w_max(ig)=0. | ||
488 | !! ale_bl_stat(ig)=0. | ||
489 | !! endif | ||
490 | !! enddo | ||
491 | 480 | susqr2pi=su*sqrt(2.*Rpi) | |
492 | reuler=exp(1.) | ||
493 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | do ig=1,ngrid |
494 |
4/4✓ Branch 0 taken 136195 times.
✓ Branch 1 taken 340925 times.
✓ Branch 2 taken 125602 times.
✓ Branch 3 taken 10593 times.
|
477600 | if ( (depth(ig).ge.10.) .and. (s_max(ig).gt.susqr2pi*reuler) ) then |
495 | 125602 | w_max(ig)=w0(ig)*(1.+sqrt(2.*log(s_max(ig)/susqr2pi)-log(2.*log(s_max(ig)/susqr2pi)))) | |
496 | 125602 | ale_bl_stat(ig)=0.5*w_max(ig)**2 | |
497 | else | ||
498 | 351518 | w_max(ig)=0. | |
499 | 351518 | ale_bl_stat(ig)=0. | |
500 | endif | ||
501 | enddo | ||
502 | |||
503 | ENDIF ! iflag_trig_bl | ||
504 | ! print *,'ENDIF iflag_trig_bl' !!jyg | ||
505 | |||
506 | !------------Closure------------------ | ||
507 | |||
508 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 480 times.
|
480 | IF (iflag_clos_bl.ge.2) THEN |
509 | |||
510 | !-----Calcul de ALP_BL_STAT | ||
511 | ✗ | do ig=1,ngrid | |
512 | ✗ | alp_bl_det(ig)=0.5*coef_m*rhobarz0(ig)*(w0(ig)**3)*fraca0(ig)*(1.-2.*fraca0(ig))/((1.-fraca0(ig))**2) | |
513 | alp_bl_fluct_m(ig)=1.5*rhobarz0(ig)*fraca0(ig)*(w_conv(ig)+coef_m*w0(ig))* & | ||
514 | ✗ | & (w0(ig)**2) | |
515 | alp_bl_fluct_tke(ig)=3.*coef_m*rhobarz0(ig)*w0(ig)*fraca0(ig)*(therm_tke_max0(ig)-env_tke_max0(ig)) & | ||
516 | ✗ | & +3.*rhobarz0(ig)*w_conv(ig)*pbl_tke_max0(ig) | |
517 | if (iflag_clos_bl.ge.2) then | ||
518 | alp_bl_conv(ig)=1.5*coef_m*rhobarz0(ig)*fraca0(ig)*(fraca0(ig)/(1.-fraca0(ig)))*w_conv(ig)* & | ||
519 | ✗ | & (w0(ig)**2) | |
520 | else | ||
521 | alp_bl_conv(ig)=0. | ||
522 | endif | ||
523 | ✗ | alp_bl_stat(ig)=alp_bl_det(ig)+alp_bl_fluct_m(ig)+alp_bl_fluct_tke(ig)+alp_bl_conv(ig) | |
524 | enddo | ||
525 | |||
526 | !-----S�curit� ALP infinie | ||
527 | ✗ | do ig=1,ngrid | |
528 | ✗ | if (fraca0(ig).gt.0.98) alp_bl_stat(ig)=2. | |
529 | enddo | ||
530 | |||
531 | ENDIF ! (iflag_clos_bl.ge.2) | ||
532 | |||
533 | !!! fin nrlmd le 10/04/2012 | ||
534 | |||
535 | ! print*,'avant calcul ale et alp' | ||
536 | !calcul de ALE et ALP pour la convection | ||
537 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | alp_bl(:)=0. |
538 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | ale_bl(:)=0. |
539 | ! print*,'ALE,ALP ,l,zw2(ig,l),ale_bl(ig),alp_bl(ig)' | ||
540 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | do l=1,nlay |
541 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | do ig=1,ngrid |
542 | 18607680 | alp_bl(ig)=max(alp_bl(ig),0.5*rhobarz(ig,l)*wth3(ig,l) ) | |
543 | 18626400 | ale_bl(ig)=max(ale_bl(ig),0.5*zw2(ig,l)**2) | |
544 | ! print*,'ALE,ALP',l,zw2(ig,l),ale_bl(ig),alp_bl(ig) | ||
545 | enddo | ||
546 | enddo | ||
547 | |||
548 | ! ale sec (max de wmax/2 sous la zone d'inhibition) dans | ||
549 | ! le cas iflag_trig_bl=3 | ||
550 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 480 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
480 | IF (iflag_trig_bl==3) ale_bl(:)=0.5*wmax_sec(:)**2 |
551 | |||
552 | !test:calcul de la ponderation des couches pour KE | ||
553 | !initialisations | ||
554 | |||
555 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | fm_tot(:)=0. |
556 |
4/4✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
✓ Branch 2 taken 18607680 times.
✓ Branch 3 taken 18720 times.
|
18626880 | wght_th(:,:)=1. |
557 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | lalim_conv(:)=lalim(:) |
558 | |||
559 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | do k=1,klev |
560 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | do ig=1,ngrid |
561 |
2/2✓ Branch 0 taken 1086781 times.
✓ Branch 1 taken 17520899 times.
|
18626400 | if (k<=lalim_conv(ig)) fm_tot(ig)=fm_tot(ig)+fm(ig,k) |
562 | enddo | ||
563 | enddo | ||
564 | |||
565 | ! assez bizarre car, si on est dans la couche d'alim et que alim_star et | ||
566 | ! plus petit que 1.e-10, on prend wght_th=1. | ||
567 |
2/2✓ Branch 0 taken 18720 times.
✓ Branch 1 taken 480 times.
|
19200 | do k=1,klev |
568 |
2/2✓ Branch 0 taken 18607680 times.
✓ Branch 1 taken 18720 times.
|
18626880 | do ig=1,ngrid |
569 |
4/4✓ Branch 0 taken 1086781 times.
✓ Branch 1 taken 17520899 times.
✓ Branch 2 taken 608736 times.
✓ Branch 3 taken 478045 times.
|
18626400 | if (k<=lalim_conv(ig).and.alim_star(ig,k)>1.e-10) then |
570 | 608736 | wght_th(ig,k)=alim_star(ig,k) | |
571 | endif | ||
572 | enddo | ||
573 | enddo | ||
574 | |||
575 | ! print*,'apres wght_th' | ||
576 | !test pour prolonger la convection | ||
577 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | do ig=1,ngrid |
578 | !v1d if ((alim_star(ig,1).lt.1.e-10).and.(therm)) then | ||
579 |
2/2✓ Branch 0 taken 228963 times.
✓ Branch 1 taken 248157 times.
|
477600 | if ((alim_star(ig,1).lt.1.e-10)) then |
580 | 228963 | lalim_conv(ig)=1 | |
581 | 228963 | wght_th(ig,1)=1. | |
582 | ! print*,'lalim_conv ok',lalim_conv(ig),wght_th(ig,1) | ||
583 | endif | ||
584 | enddo | ||
585 | |||
586 | !------------------------------------------------------------------------ | ||
587 | ! Modif CR/FH 20110310 : alp integree sur la verticale. | ||
588 | ! Integrale verticale de ALP. | ||
589 | ! wth3 etant aux niveaux inter-couches, on utilise d play comme masse des | ||
590 | ! couches | ||
591 | !------------------------------------------------------------------------ | ||
592 | |||
593 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | alp_int(:)=0. |
594 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | dp_int(:)=0. |
595 |
2/2✓ Branch 0 taken 18240 times.
✓ Branch 1 taken 480 times.
|
18720 | do l=2,nlay |
596 |
2/2✓ Branch 0 taken 18130560 times.
✓ Branch 1 taken 18240 times.
|
18149280 | do ig=1,ngrid |
597 |
2/2✓ Branch 0 taken 1306490 times.
✓ Branch 1 taken 16824070 times.
|
18148800 | if(l.LE.lmax(ig)) THEN |
598 | 1306490 | zdp=pplay(ig,l-1)-pplay(ig,l) | |
599 | 1306490 | alp_int(ig)=alp_int(ig)+0.5*rhobarz(ig,l)*wth3(ig,l)*zdp | |
600 | 1306490 | dp_int(ig)=dp_int(ig)+zdp | |
601 | endif | ||
602 | enddo | ||
603 | enddo | ||
604 | |||
605 |
1/2✓ Branch 0 taken 480 times.
✗ Branch 1 not taken.
|
480 | if (iflag_coupl>=3 .and. iflag_coupl<=5) then |
606 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | do ig=1,ngrid |
607 | !valeur integree de alp_bl * 0.5: | ||
608 |
2/2✓ Branch 0 taken 248430 times.
✓ Branch 1 taken 228690 times.
|
477600 | if (dp_int(ig)>0.) then |
609 | 248430 | alp_bl(ig)=alp_int(ig)/dp_int(ig) | |
610 | endif | ||
611 | enddo! | ||
612 | endif | ||
613 | |||
614 | |||
615 | ! Facteur multiplicatif sur alp_bl | ||
616 |
2/2✓ Branch 0 taken 477120 times.
✓ Branch 1 taken 480 times.
|
477600 | alp_bl(:)=alp_bl_k*alp_bl(:) |
617 | |||
618 | !------------------------------------------------------------------------ | ||
619 | |||
620 | |||
621 | |||
622 | 480 | return | |
623 | end | ||
624 |