/[lmdze]/trunk/phylmd/hbtm.f
ViewVC logotype

Contents of /trunk/phylmd/hbtm.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 254 - (show annotations)
Mon Feb 5 10:39:38 2018 UTC (6 years, 3 months ago) by guez
File size: 13346 byte(s)
Move Sources/* to root directory.
1 module HBTM_m
2
3 IMPLICIT none
4
5 contains
6
7 SUBROUTINE HBTM(paprs, pplay, t2m, q2m, ustar, flux_t, flux_q, u, v, t, q, &
8 pblh, cape, EauLiq, ctei, pblT, therm, plcl)
9
10 ! D'apr\'es Holstag et Boville et Troen et Mahrt
11 ! JAS 47 BLM
12
13 ! Algorithme th\'ese Anne Mathieu. Crit\'ere d'entra\^inement
14 ! Peter Duynkerke (JAS 50). Written by: Anne MATHIEU and Alain
15 ! LAHELLEC, 22nd November 1999.
16
17 ! Modifications : d\'ecembre 99 passage th \`a niveau plus bas. Voir fixer
18 ! la prise du th \`a z/Lambda = -.2 (max Ray)
19 ! Autre algorithme : entra\^inement ~ Theta + v =constante
20 ! mais comment ? The ?
21 ! On peut fixer q \`a 0.7 qsat (cf. non adiabatique) d'où T2 et The2.
22 ! Voir aussi KE pblh = niveau The_e ou l = env.
23
24 ! Adaptation \`a LMDZ version coupl\'ee. Pour le moment on fait
25 ! passer en argument les grandeurs de surface : flux, t, q2m. On
26 ! va utiliser syst\'ematiquement les grandeurs \`a 2 m mais on
27 ! garde la possibilit\'e de changer si besoin (jusqu'\`a pr\'esent
28 ! la forme de HB avec le premier niveau mod\`ele \'etait
29 ! conserv\'ee).
30
31 USE dimphy, ONLY: klev, klon
32 USE suphec_m, ONLY: rcpd, rd, retv, rg, rkappa, rtt
33 USE yoethf_m, ONLY: r2es, rvtmp2
34 USE fcttre, ONLY: foeew
35
36 ! Arguments:
37
38 ! pression a inter-couche (Pa)
39 REAL, intent(in):: paprs(klon, klev+1)
40 ! pression au milieu de couche (Pa)
41 REAL, intent(in):: pplay(klon, klev)
42 REAL, intent(in):: t2m(klon) ! temperature a 2 m
43 ! q a 2 et 10m
44 REAL, intent(in):: q2m(klon)
45 REAL, intent(in):: ustar(:) ! (knon)
46 REAL, intent(in):: flux_t(:), flux_q(:) ! (knon) flux à la surface
47
48 REAL, intent(in):: u(klon, klev) ! vitesse U (m/s)
49 REAL, intent(in):: v(klon, klev) ! vitesse V (m/s)
50
51 ! temperature (K)
52 REAL, intent(in):: t(klon, klev)
53 ! vapeur d'eau (kg/kg)
54 REAL, intent(in):: q(klon, klev)
55
56 REAL, intent(out):: pblh(:) ! (knon)
57 ! Cape du thermique
58 REAL Cape(klon)
59 ! Eau liqu integr du thermique
60 REAL EauLiq(klon)
61 ! Critere d'instab d'entrainmt des nuages de
62 REAL ctei(klon)
63 REAL pblT(klon)
64 ! thermal virtual temperature excess
65 REAL therm(klon)
66 REAL plcl(klon)
67
68 ! Local:
69
70 INTEGER knon ! nombre de points a calculer
71 INTEGER isommet
72 ! limite max sommet pbl
73 PARAMETER (isommet=klev)
74 REAL vk
75 ! Von Karman => passer a .41 ! cf U.Olgstrom
76 PARAMETER (vk=0.35)
77 REAL ricr
78 PARAMETER (ricr=0.4)
79 ! a
80 REAL onet
81 PARAMETER (onet=1.0/3.0)
82 REAL zkmin
83 PARAMETER (zkmin=0.01)
84 REAL betam
85 ! pour Phim / h dans la S.L stable
86 PARAMETER (betam=15.0)
87 ! z/OBL<>1
88 REAL sffrac
89 ! S.L. = z/h < .1
90 PARAMETER (sffrac=0.1)
91 REAL binm
92 PARAMETER (binm=betam*sffrac)
93
94 REAL q_star, t_star
95 ! Lambert correlations T' q' avec T* q*
96 REAL b1, b2, b212, b2sr
97 PARAMETER (b1=70., b2=20.)
98
99 REAL z(klon, klev)
100
101 REAL zref
102 ! Niveau de ref a 2m peut eventuellement
103 PARAMETER (zref=2.)
104 ! etre choisi a 10m
105
106 INTEGER i, k
107 REAL zxt
108 ! surface kinematic heat flux [mK/s]
109 REAL khfs(klon)
110 ! sfc kinematic constituent flux [m/s]
111 REAL kqfs(klon)
112 ! surface virtual heat flux
113 REAL heatv(klon)
114 ! bulk Richardon no. mais en Theta_v
115 REAL rhino(klon, klev)
116 ! pts w/unstbl pbl (positive virtual ht flx)
117 LOGICAL unstbl(klon)
118 LOGICAL check(klon) ! Richardson number > critical
119 ! flag de prolongerment cape pour pt Omega
120 LOGICAL omegafl(klon)
121
122 ! Monin-Obukhov lengh
123 REAL obklen(klon)
124
125 REAL zdu2
126 ! Algorithme thermique
127 REAL s(klon, klev) ! [P/Po]^Kappa milieux couches
128 ! total water of thermal
129 REAL qT_th(klon)
130 ! T thermique niveau precedent
131 REAL qsatbef(klon)
132 ! le thermique est sature
133 LOGICAL Zsat(klon)
134 REAL zthvd, zthvu, qqsat
135 REAL t2
136
137 ! inverse phi function for momentum
138 REAL phiminv(klon)
139 ! turbulent velocity scale for momentum
140 REAL wm(klon)
141 ! current level height + one level up
142 REAL zp(klon)
143 REAL zcor
144
145 REAL pblmin
146
147 !-----------------------------------------------------------------
148
149 knon = size(pblh)
150
151 ! initialisations
152 q_star = 0
153 t_star = 0
154
155 b212=sqrt(b1*b2)
156 b2sr=sqrt(b2)
157
158 ! Calculer les hauteurs de chaque couche
159 ! (geopotentielle Int_dp/ro = Int_[Rd.T.dp/p] z = geop/g)
160 ! pourquoi ne pas utiliser Phi/RG ?
161 DO i = 1, knon
162 z(i, 1) = RD * t(i, 1) / (0.5*(paprs(i, 1)+pplay(i, 1))) &
163 * (paprs(i, 1)-pplay(i, 1)) / RG
164 s(i, 1) = (pplay(i, 1)/paprs(i, 1))**RKappa
165 ENDDO
166 ! s(k) = [pplay(k)/ps]^kappa
167 ! + + + + + + + + + pplay <-> s(k) t dp=pplay(k-1)-pplay(k)
168 ! ----------------- paprs <-> sig(k)
169 ! + + + + + + + + + pplay <-> s(k-1)
170 ! + + + + + + + + + pplay <-> s(1) t dp=paprs-pplay z(1)
171 ! ----------------- paprs <-> sig(1)
172
173 DO k = 2, klev
174 DO i = 1, knon
175 z(i, k) = z(i, k-1) &
176 + RD * 0.5*(t(i, k-1)+t(i, k)) / paprs(i, k) &
177 * (pplay(i, k-1)-pplay(i, k)) / RG
178 s(i, k) = (pplay(i, k) / paprs(i, 1))**RKappa
179 ENDDO
180 ENDDO
181
182 ! Determination des grandeurs de surface
183 DO i = 1, knon
184 ! Niveau de ref choisi a 2m
185 zxt = t2m(i)
186
187 ! convention >0 vers le bas ds lmdz
188 khfs(i) = - flux_t(i)*zxt*Rd / (RCPD*paprs(i, 1))
189 kqfs(i) = - flux_q(i)*zxt*Rd / paprs(i, 1)
190 ! verifier que khfs et kqfs sont bien de la forme w'l'
191 heatv(i) = khfs(i) + 0.608*zxt*kqfs(i)
192 ! a comparer aussi aux sorties de clqh : flux_T/RoCp et flux_q/RoLv
193 ! Theta et qT du thermique sans exces (interpolin vers surf)
194 ! chgt de niveau du thermique (jeudi 30/12/1999)
195 ! (interpolation lineaire avant integration phi_h)
196 qT_th(i) = q2m(i)
197 ENDDO
198
199 DO i = 1, knon
200 ! Global Richardson
201 rhino(i, 1) = 0.0
202 check(i) = .TRUE.
203 ! on initialise pblh a l'altitude du 1er niv
204 pblh(i) = z(i, 1)
205 plcl(i) = 6000.
206 ! Lambda = -u*^3 / (alpha.g.kvon.<w'Theta'v>
207 obklen(i) = -t(i, 1)*ustar(i)**3/(RG*vk*heatv(i))
208 ENDDO
209
210 ! PBL height calculation: Search for level of pbl. Scan upward
211 ! until the Richardson number between the first level and the
212 ! current level exceeds the "critical" value. (bonne idee Nu de
213 ! separer le Ric et l'exces de temp du thermique)
214 DO k = 2, isommet
215 DO i = 1, knon
216 IF (check(i)) THEN
217 ! pourquoi / niveau 1 (au lieu du sol) et le terme en u*^2 ?
218 zdu2 = u(i, k)**2+v(i, k)**2
219 zdu2 = max(zdu2, 1.0e-20)
220 ! Theta_v environnement
221 zthvd=t(i, k)/s(i, k)*(1.+RETV*q(i, k))
222
223 ! therm Theta_v sans exces (avec hypothese fausse de H&B, sinon,
224 ! passer par Theta_e et virpot)
225 zthvu = T2m(i)*(1.+RETV*qT_th(i))
226 ! Le Ri par Theta_v
227 ! On a nveau de ref a 2m ???
228 rhino(i, k) = (z(i, k)-zref)*RG*(zthvd-zthvu) &
229 /(zdu2*0.5*(zthvd+zthvu))
230
231 IF (rhino(i, k).GE.ricr) THEN
232 pblh(i) = z(i, k-1) + (z(i, k-1)-z(i, k)) * &
233 (ricr-rhino(i, k-1))/(rhino(i, k-1)-rhino(i, k))
234 ! test04
235 pblh(i) = pblh(i) + 100.
236 pblT(i) = t(i, k-1) + (t(i, k)-t(i, k-1)) * &
237 (pblh(i)-z(i, k-1))/(z(i, k)-z(i, k-1))
238 check(i) = .FALSE.
239 ENDIF
240 ENDIF
241 ENDDO
242 ENDDO
243
244 ! Set pbl height to maximum value where computation exceeds number of
245 ! layers allowed
246 DO i = 1, knon
247 if (check(i)) pblh(i) = z(i, isommet)
248 ENDDO
249
250 ! Improve estimate of pbl height for the unstable points.
251 ! Find unstable points (sensible heat flux is upward):
252 DO i = 1, knon
253 IF (heatv(i) > 0.) THEN
254 unstbl(i) = .TRUE.
255 check(i) = .TRUE.
256 ELSE
257 unstbl(i) = .FALSE.
258 check(i) = .FALSE.
259 ENDIF
260 ENDDO
261
262 ! For the unstable case, compute velocity scale and the
263 ! convective temperature excess:
264 DO i = 1, knon
265 IF (check(i)) THEN
266 phiminv(i) = (1.-binm*pblh(i)/obklen(i))**onet
267
268 ! CALCUL DE wm
269 ! Ici on considerera que l'on est dans la couche de surf jusqu'a 100
270 ! On prend svt couche de surface=0.1*h mais on ne connait pas h
271 ! Dans la couche de surface
272 wm(i)= ustar(i)*phiminv(i)
273
274 ! forme Mathieu :
275 q_star = kqfs(i)/wm(i)
276 t_star = khfs(i)/wm(i)
277
278 therm(i) = sqrt( b1*(1.+2.*RETV*qT_th(i))*t_star**2 &
279 + (RETV*T2m(i))**2*b2*q_star**2 &
280 + max(0., 2.*RETV*T2m(i)*b212*q_star*t_star))
281
282 ! Theta et qT du thermique (forme H&B) avec exces
283 ! (attention, on ajoute therm(i) qui est virtuelle ...)
284 ! pourquoi pas sqrt(b1)*t_star ?
285 qT_th(i) = qT_th(i) + b2sr*q_star
286 ! new on differre le calcul de Theta_e
287 rhino(i, 1) = 0.
288 ENDIF
289 ENDDO
290
291 ! Improve pblh estimate for unstable conditions using the
292 ! convective temperature excess :
293 DO k = 2, isommet
294 DO i = 1, knon
295 IF (check(i)) THEN
296 zdu2 = u(i, k)**2 + v(i, k)**2
297 zdu2 = max(zdu2, 1e-20)
298 ! Theta_v environnement
299 zthvd=t(i, k)/s(i, k)*(1.+RETV*q(i, k))
300
301 ! et therm Theta_v (avec hypothese de constance de H&B,
302 zthvu = T2m(i)*(1.+RETV*qT_th(i)) + therm(i)
303
304 ! Le Ri par Theta_v
305 ! Niveau de ref 2m
306 rhino(i, k) = (z(i, k)-zref)*RG*(zthvd-zthvu) &
307 /(zdu2*0.5*(zthvd+zthvu))
308
309 IF (rhino(i, k).GE.ricr) THEN
310 pblh(i) = z(i, k-1) + (z(i, k-1)-z(i, k)) * &
311 (ricr-rhino(i, k-1))/(rhino(i, k-1)-rhino(i, k))
312 ! test04
313 pblh(i) = pblh(i) + 100.
314 pblT(i) = t(i, k-1) + (t(i, k)-t(i, k-1)) * &
315 (pblh(i)-z(i, k-1))/(z(i, k)-z(i, k-1))
316 check(i) = .FALSE.
317 ENDIF
318 ENDIF
319 ENDDO
320 ENDDO
321
322 ! Set pbl height to maximum value where computation exceeds number of
323 ! layers allowed
324 DO i = 1, knon
325 if (check(i)) pblh(i) = z(i, isommet)
326 ENDDO
327
328 ! PBL height must be greater than some minimum mechanical mixing depth
329 ! Several investigators have proposed minimum mechanical mixing depth
330 ! relationships as a function of the local friction velocity, u*. We
331 ! make use of a linear relationship of the form h = c u* where c=700.
332 ! The scaling arguments that give rise to this relationship most often
333 ! represent the coefficient c as some constant over the local coriolis
334 ! parameter. Here we make use of the experimental results of Koracin
335 ! and Berkowicz (1988) [BLM, Vol 43] for wich they recommend 0.07/f
336 ! where f was evaluated at 39.5 N and 52 N. Thus we use a typical mid
337 ! latitude value for f so that c = 0.07/f = 700.
338 DO i = 1, knon
339 pblmin = 700. * ustar(i)
340 pblh(i) = MAX(pblh(i), pblmin)
341 ! par exemple :
342 pblT(i) = t(i, 2) + (t(i, 3)-t(i, 2)) * &
343 (pblh(i)-z(i, 2))/(z(i, 3)-z(i, 2))
344 ENDDO
345
346 ! pblh is now available; do preparation for diffusivity calculation:
347 DO i = 1, knon
348 check(i) = .TRUE.
349 Zsat(i) = .FALSE.
350 ! omegafl utilise pour prolongement CAPE
351 omegafl(i) = .FALSE.
352 Cape(i) = 0.
353 EauLiq(i) = 0.
354 CTEI(i) = 0.
355
356 ! Do additional preparation for unstable cases only, set temperature
357 ! and moisture perturbations depending on stability.
358 ! Remarque : les formule sont prises dans leur forme CS
359 IF (unstbl(i)) THEN
360 ! Niveau de ref du thermique
361 zxt=(T2m(i)-zref*0.5*RG/RCPD/(1.+RVTMP2*qT_th(i))) &
362 *(1.+RETV*qT_th(i))
363 phiminv(i) = (1. - binm*pblh(i)/obklen(i))**onet
364 wm(i) = ustar(i)*phiminv(i)
365 ENDIF
366 ENDDO
367
368 ! Main level loop to compute the diffusivities and
369 ! counter-gradient terms:
370 loop_level: DO k = 2, isommet
371 ! Find levels within boundary layer:
372 DO i = 1, knon
373 zp(i) = z(i, k)
374 IF (zkmin == 0. .AND. zp(i) > pblh(i)) zp(i) = pblh(i)
375 ENDDO
376
377 ! For all layers, compute integral info and CTEI
378 DO i = 1, knon
379 if (check(i) .or. omegafl(i)) then
380 if (.not. Zsat(i)) then
381 T2 = T2m(i) * s(i, k)
382 ! thermodyn functions
383 qqsat= r2es * FOEEW(T2, RTT >= T2) / pplay(i, k)
384 qqsat=MIN(0.5, qqsat)
385 zcor=1./(1.-retv*qqsat)
386 qqsat=qqsat*zcor
387
388 if (qqsat < qT_th(i)) then
389 ! on calcule lcl
390 if (k == 2) then
391 plcl(i) = z(i, k)
392 else
393 plcl(i) = z(i, k-1) + (z(i, k-1)-z(i, k)) &
394 * (qT_th(i)-qsatbef(i)) / (qsatbef(i)-qqsat)
395 endif
396 Zsat(i) = .true.
397 endif
398 endif
399 qsatbef(i) = qqsat
400 ! cette ligne a deja ete faite normalement ?
401 endif
402 ENDDO
403 end DO loop_level
404
405 END SUBROUTINE HBTM
406
407 end module HBTM_m

Properties

Name Value
svn:executable

  ViewVC Help
Powered by ViewVC 1.1.21