1 | SUBROUTINE HBTM(knon, paprs, pplay, |
---|
2 | . t2m,t10m,q2m,q10m,ustar, |
---|
3 | . flux_t,flux_q,u,v,t,q, |
---|
4 | . pblh,cape,EauLiq,ctei,pblT, |
---|
5 | . therm,trmb1,trmb2,trmb3,plcl) |
---|
6 | IMPLICIT none |
---|
7 | |
---|
8 | c*************************************************************** |
---|
9 | c* * |
---|
10 | c* HBTM2 D'apres Holstag&Boville et Troen&Mahrt * |
---|
11 | c* JAS 47 BLM * |
---|
12 | c* Algorithme These Anne Mathieu * |
---|
13 | c* Critere d'Entrainement Peter Duynkerke (JAS 50) * |
---|
14 | c* written by : Anne MATHIEU & Alain LAHELLEC, 22/11/99 * |
---|
15 | c* features : implem. exces Mathieu * |
---|
16 | c*************************************************************** |
---|
17 | c* mods : decembre 99 passage th a niveau plus bas. voir fixer * |
---|
18 | c* la prise du th a z/Lambda = -.2 (max Ray) * |
---|
19 | c* Autre algo : entrainement ~ Theta+v =cste mais comment=>The?* |
---|
20 | c* on peut fixer q a .7qsat(cf non adiab)=>T2 et The2 * |
---|
21 | c* voir aussi //KE pblh = niveau The_e ou l = env. * |
---|
22 | c*************************************************************** |
---|
23 | c* fin therm a la HBTM passage a forme Mathieu 12/09/2001 * |
---|
24 | c*************************************************************** |
---|
25 | c* |
---|
26 | c |
---|
27 | c |
---|
28 | cAM Fev 2003 |
---|
29 | c Adaptation a LMDZ version couplee |
---|
30 | c |
---|
31 | c Pour le moment on fait passer en argument les grdeurs de surface : |
---|
32 | c flux, t,q2m, t,q10m, on va utiliser systematiquement les grdeurs a 2m ms |
---|
33 | c on garde la possibilite de changer si besoin est (jusqu'a present la |
---|
34 | c forme de HB avec le 1er niveau modele etait conservee) |
---|
35 | c |
---|
36 | c |
---|
37 | c |
---|
38 | c |
---|
39 | c |
---|
40 | #include "dimensions.h" |
---|
41 | #include "dimphy.h" |
---|
42 | #include "YOMCST.h" |
---|
43 | REAL RLvCp, REPS |
---|
44 | c Arguments: |
---|
45 | c |
---|
46 | INTEGER knon ! nombre de points a calculer |
---|
47 | cAM |
---|
48 | REAL t2m(klon), t10m(klon) ! temperature a 2 et 10m |
---|
49 | REAL q2m(klon), q10m(klon) ! q a 2 et 10m |
---|
50 | REAL ustar(klon) |
---|
51 | REAL paprs(klon,klev+1) ! pression a inter-couche (Pa) |
---|
52 | REAL pplay(klon,klev) ! pression au milieu de couche (Pa) |
---|
53 | REAL flux_t(klon,klev), flux_q(klon,klev) ! Flux |
---|
54 | REAL u(klon,klev) ! vitesse U (m/s) |
---|
55 | REAL v(klon,klev) ! vitesse V (m/s) |
---|
56 | REAL t(klon,klev) ! temperature (K) |
---|
57 | REAL q(klon,klev) ! vapeur d'eau (kg/kg) |
---|
58 | cAM REAL cd_h(klon) ! coefficient de friction au sol pour chaleur |
---|
59 | cAM REAL cd_m(klon) ! coefficient de friction au sol pour vitesse |
---|
60 | c |
---|
61 | INTEGER isommet |
---|
62 | PARAMETER (isommet=klev) ! limite max sommet pbl |
---|
63 | REAL vk |
---|
64 | PARAMETER (vk=0.35) ! Von Karman => passer a .41 ! cf U.Olgstrom |
---|
65 | REAL ricr |
---|
66 | PARAMETER (ricr=0.4) |
---|
67 | REAL fak |
---|
68 | PARAMETER (fak=8.5) ! b calcul du Prandtl et de dTetas |
---|
69 | REAL fakn |
---|
70 | PARAMETER (fakn=7.2) ! a |
---|
71 | REAL onet |
---|
72 | PARAMETER (onet=1.0/3.0) |
---|
73 | REAL t_coup |
---|
74 | PARAMETER(t_coup=273.15) |
---|
75 | REAL zkmin |
---|
76 | PARAMETER (zkmin=0.01) |
---|
77 | REAL betam |
---|
78 | PARAMETER (betam=15.0) ! pour Phim / h dans la S.L stable |
---|
79 | REAL betah |
---|
80 | PARAMETER (betah=15.0) |
---|
81 | REAL betas |
---|
82 | PARAMETER (betas=5.0) ! Phit dans la S.L. stable (mais 2 formes / z/OBL<>1 |
---|
83 | REAL sffrac |
---|
84 | PARAMETER (sffrac=0.1) ! S.L. = z/h < .1 |
---|
85 | REAL binm |
---|
86 | PARAMETER (binm=betam*sffrac) |
---|
87 | REAL binh |
---|
88 | PARAMETER (binh=betah*sffrac) |
---|
89 | REAL ccon |
---|
90 | PARAMETER (ccon=fak*sffrac*vk) |
---|
91 | c |
---|
92 | REAL q_star,t_star |
---|
93 | REAL b1,b2,b212,b2sr ! Lambert correlations T' q' avec T* q* |
---|
94 | PARAMETER (b1=70.,b2=20.) |
---|
95 | c |
---|
96 | REAL z(klon,klev) |
---|
97 | cAM REAL pcfm(klon,klev), pcfh(klon,klev) |
---|
98 | cAM |
---|
99 | REAL zref |
---|
100 | PARAMETER (zref=2.) ! Niveau de ref a 2m peut eventuellement |
---|
101 | c etre choisi a 10m |
---|
102 | cMA |
---|
103 | c |
---|
104 | INTEGER i, k, j |
---|
105 | REAL zxt |
---|
106 | cAM REAL zxt, zxq, zxu, zxv, zxmod, taux, tauy |
---|
107 | cAM REAL zx_alf1, zx_alf2 ! parametres pour extrapolation |
---|
108 | REAL khfs(klon) ! surface kinematic heat flux [mK/s] |
---|
109 | REAL kqfs(klon) ! sfc kinematic constituent flux [m/s] |
---|
110 | REAL heatv(klon) ! surface virtual heat flux |
---|
111 | REAL rhino(klon,klev) ! bulk Richardon no. mais en Theta_v |
---|
112 | LOGICAL unstbl(klon) ! pts w/unstbl pbl (positive virtual ht flx) |
---|
113 | LOGICAL stblev(klon) ! stable pbl with levels within pbl |
---|
114 | LOGICAL unslev(klon) ! unstbl pbl with levels within pbl |
---|
115 | LOGICAL unssrf(klon) ! unstb pbl w/lvls within srf pbl lyr |
---|
116 | LOGICAL unsout(klon) ! unstb pbl w/lvls in outer pbl lyr |
---|
117 | LOGICAL check(klon) ! True=>chk if Richardson no.>critcal |
---|
118 | LOGICAL omegafl(klon) ! flag de prolongerment cape pour pt Omega |
---|
119 | REAL pblh(klon) |
---|
120 | REAL pblT(klon) |
---|
121 | REAL plcl(klon) |
---|
122 | cAM REAL cgh(klon,2:klev) ! counter-gradient term for heat [K/m] |
---|
123 | cAM REAL cgq(klon,2:klev) ! counter-gradient term for constituents |
---|
124 | cAM REAL cgs(klon,2:klev) ! counter-gradient star (cg/flux) |
---|
125 | REAL obklen(klon) ! Monin-Obukhov lengh |
---|
126 | cAM REAL ztvd, ztvu, |
---|
127 | REAL zdu2 |
---|
128 | REAL therm(klon) ! thermal virtual temperature excess |
---|
129 | REAL trmb1(klon),trmb2(klon),trmb3(klon) |
---|
130 | C Algorithme thermique |
---|
131 | REAL s(klon,klev) ! [P/Po]^Kappa milieux couches |
---|
132 | REAL Th_th(klon) ! potential temperature of thermal |
---|
133 | REAL The_th(klon) ! equivalent potential temperature of thermal |
---|
134 | REAL qT_th(klon) ! total water of thermal |
---|
135 | REAL Tbef(klon) ! T thermique niveau precedent |
---|
136 | REAL qsatbef(klon) |
---|
137 | LOGICAL Zsat(klon) ! le thermique est sature |
---|
138 | REAL Cape(klon) ! Cape du thermique |
---|
139 | REAL Kape(klon) ! Cape locale |
---|
140 | REAL EauLiq(klon) ! Eau liqu integr du thermique |
---|
141 | REAL ctei(klon) ! Critere d'instab d'entrainmt des nuages de CL |
---|
142 | REAL the1,the2,aa,bb,zthvd,zthvu,xintpos,qqsat |
---|
143 | cIM 091204 BEG |
---|
144 | REAL a1,a2,a3 |
---|
145 | cIM 091204 END |
---|
146 | REAL xhis,rnum,denom,th1,th2,thv1,thv2,ql2 |
---|
147 | REAL dqsat_dt,qsat2,qT1,q2,t1,t2,xnull,delt_the |
---|
148 | REAL delt_qt,delt_2,quadsat,spblh,reduc |
---|
149 | c |
---|
150 | REAL phiminv(klon) ! inverse phi function for momentum |
---|
151 | REAL phihinv(klon) ! inverse phi function for heat |
---|
152 | REAL wm(klon) ! turbulent velocity scale for momentum |
---|
153 | REAL fak1(klon) ! k*ustar*pblh |
---|
154 | REAL fak2(klon) ! k*wm*pblh |
---|
155 | REAL fak3(klon) ! fakn*wstr/wm |
---|
156 | REAL pblk(klon) ! level eddy diffusivity for momentum |
---|
157 | REAL pr(klon) ! Prandtl number for eddy diffusivities |
---|
158 | REAL zl(klon) ! zmzp / Obukhov length |
---|
159 | REAL zh(klon) ! zmzp / pblh |
---|
160 | REAL zzh(klon) ! (1-(zmzp/pblh))**2 |
---|
161 | REAL wstr(klon) ! w*, convective velocity scale |
---|
162 | REAL zm(klon) ! current level height |
---|
163 | REAL zp(klon) ! current level height + one level up |
---|
164 | REAL zcor, zdelta, zcvm5 |
---|
165 | cAM REAL zxqs |
---|
166 | REAL fac, pblmin, zmzp, term |
---|
167 | c |
---|
168 | #include "YOETHF.h" |
---|
169 | #include "FCTTRE.h" |
---|
170 | |
---|
171 | |
---|
172 | |
---|
173 | ! initialisations (Anne) |
---|
174 | th_th(:) = 0. |
---|
175 | q_star = 0 |
---|
176 | t_star = 0 |
---|
177 | |
---|
178 | |
---|
179 | b212=sqrt(b1*b2) |
---|
180 | b2sr=sqrt(b2) |
---|
181 | c |
---|
182 | C ============================================================ |
---|
183 | C Fonctions thermo implicites |
---|
184 | C ============================================================ |
---|
185 | c +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
186 | c Tetens : pression partielle de vap d'eau e_sat(T) |
---|
187 | c ================================================= |
---|
188 | C++ e_sat(T) = r2*exp( r3*(T-Tf)/(T-r4) ) id a r2*FOEWE |
---|
189 | C++ avec : |
---|
190 | C++ Tf = 273.16 K (Temp de fusion de la glace) |
---|
191 | C++ r2 = 611.14 Pa |
---|
192 | C++ r3 = 17.269 (liquide) 21.875 (solide) adim |
---|
193 | C++ r4 = 35.86 7.66 Kelvin |
---|
194 | C++ q_sat = eps*e_sat/(p-(1-eps)*e_sat) |
---|
195 | C++ derivée : |
---|
196 | C++ ========= |
---|
197 | C++ r3*(Tf-r4)*q_sat(T,p) |
---|
198 | C++ d_qsat_dT = -------------------------------- |
---|
199 | C++ (T-r4)^2*( 1-(1-eps)*e_sat(T)/p ) |
---|
200 | c++ pour zcvm5=Lv, c'est FOEDE |
---|
201 | c++ Rq :(1.-REPS)*esarg/Parg id a RETV*Qsat |
---|
202 | C ------------------------------------------------------------------ |
---|
203 | c |
---|
204 | c Initialisation |
---|
205 | RLvCp = RLVTT/RCPD |
---|
206 | REPS = RD/RV |
---|
207 | |
---|
208 | c |
---|
209 | c DO i = 1, klon |
---|
210 | c pcfh(i,1) = cd_h(i) |
---|
211 | c pcfm(i,1) = cd_m(i) |
---|
212 | c ENDDO |
---|
213 | c DO k = 2, klev |
---|
214 | c DO i = 1, klon |
---|
215 | c pcfh(i,k) = zkmin |
---|
216 | c pcfm(i,k) = zkmin |
---|
217 | c cgs(i,k) = 0.0 |
---|
218 | c cgh(i,k) = 0.0 |
---|
219 | c cgq(i,k) = 0.0 |
---|
220 | c ENDDO |
---|
221 | c ENDDO |
---|
222 | c |
---|
223 | c Calculer les hauteurs de chaque couche |
---|
224 | c (geopotentielle Int_dp/ro = Int_[Rd.T.dp/p] z = geop/g) |
---|
225 | c pourquoi ne pas utiliser Phi/RG ? |
---|
226 | DO i = 1, knon |
---|
227 | z(i,1) = RD * t(i,1) / (0.5*(paprs(i,1)+pplay(i,1))) |
---|
228 | . * (paprs(i,1)-pplay(i,1)) / RG |
---|
229 | s(i,1) = (pplay(i,1)/paprs(i,1))**RKappa |
---|
230 | ENDDO |
---|
231 | c s(k) = [pplay(k)/ps]^kappa |
---|
232 | c + + + + + + + + + pplay <-> s(k) t dp=pplay(k-1)-pplay(k) |
---|
233 | c |
---|
234 | c ----------------- paprs <-> sig(k) |
---|
235 | c |
---|
236 | c + + + + + + + + + pplay <-> s(k-1) |
---|
237 | c |
---|
238 | c |
---|
239 | c + + + + + + + + + pplay <-> s(1) t dp=paprs-pplay z(1) |
---|
240 | c |
---|
241 | c ----------------- paprs <-> sig(1) |
---|
242 | c |
---|
243 | DO k = 2, klev |
---|
244 | DO i = 1, knon |
---|
245 | z(i,k) = z(i,k-1) |
---|
246 | . + RD * 0.5*(t(i,k-1)+t(i,k)) / paprs(i,k) |
---|
247 | . * (pplay(i,k-1)-pplay(i,k)) / RG |
---|
248 | s(i,k) = (pplay(i,k)/paprs(i,1))**RKappa |
---|
249 | ENDDO |
---|
250 | ENDDO |
---|
251 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
252 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
253 | c +++ Determination des grandeurs de surface +++++++++++++++++++++ |
---|
254 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
255 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
256 | DO i = 1, knon |
---|
257 | cAM IF (thermcep) THEN |
---|
258 | cAM zdelta=MAX(0.,SIGN(1.,RTT-tsol(i))) |
---|
259 | c zcvm5 = R5LES*RLVTT*(1.-zdelta) + R5IES*RLSTT*zdelta |
---|
260 | c zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q(i,1)) |
---|
261 | cAM zxqs= r2es * FOEEW(tsol(i),zdelta)/paprs(i,1) |
---|
262 | cAM zxqs=MIN(0.5,zxqs) |
---|
263 | cAM zcor=1./(1.-retv*zxqs) |
---|
264 | cAM zxqs=zxqs*zcor |
---|
265 | cAM ELSE |
---|
266 | cAM IF (tsol(i).LT.t_coup) THEN |
---|
267 | cAM zxqs = qsats(tsol(i)) / paprs(i,1) |
---|
268 | cAM ELSE |
---|
269 | cAM zxqs = qsatl(tsol(i)) / paprs(i,1) |
---|
270 | cAM ENDIF |
---|
271 | cAM ENDIF |
---|
272 | c niveau de reference bulk; mais ici, c,a pourrait etre le niveau de ref du thermique |
---|
273 | cAM zx_alf1 = 1.0 |
---|
274 | cAM zx_alf2 = 1.0 - zx_alf1 |
---|
275 | cAM zxt = (t(i,1)+z(i,1)*RG/RCPD/(1.+RVTMP2*q(i,1))) |
---|
276 | cAM . *(1.+RETV*q(i,1))*zx_alf1 |
---|
277 | cAM . + (t(i,2)+z(i,2)*RG/RCPD/(1.+RVTMP2*q(i,2))) |
---|
278 | cAM . *(1.+RETV*q(i,2))*zx_alf2 |
---|
279 | cAM zxu = u(i,1)*zx_alf1+u(i,2)*zx_alf2 |
---|
280 | cAM zxv = v(i,1)*zx_alf1+v(i,2)*zx_alf2 |
---|
281 | cAM zxq = q(i,1)*zx_alf1+q(i,2)*zx_alf2 |
---|
282 | cAM |
---|
283 | cAMAM zxu = u10m(i) |
---|
284 | cAMAM zxv = v10m(i) |
---|
285 | cAMAM zxmod = 1.0+SQRT(zxu**2+zxv**2) |
---|
286 | cAM Niveau de ref choisi a 2m |
---|
287 | zxt = t2m(i) |
---|
288 | |
---|
289 | c *************************************************** |
---|
290 | c attention, il doit s'agir de <w'theta'> |
---|
291 | c ;Calcul de tcls virtuel et de w'theta'virtuel |
---|
292 | c ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;; |
---|
293 | c tcls=tcls*(1+.608*qcls) |
---|
294 | c |
---|
295 | c ;Pour avoir w'theta', |
---|
296 | c ; il faut diviser par ro.Cp |
---|
297 | c Cp=Cpd*(1+0.84*qcls) |
---|
298 | c fcs=fcs/(ro_surf*Cp) |
---|
299 | c ;On transforme w'theta' en w'thetav' |
---|
300 | c Lv=(2.501-0.00237*(tcls-273.15))*1.E6 |
---|
301 | c xle=xle/(ro_surf*Lv) |
---|
302 | c fcsv=fcs+.608*xle*tcls |
---|
303 | c *************************************************** |
---|
304 | cAM khfs(i) = (tsol(i)*(1.+RETV*q(i,1))-zxt) *zxmod*cd_h(i) |
---|
305 | cAM kqfs(i) = (zxqs-zxq) *zxmod*cd_h(i) * beta(i) |
---|
306 | cAM |
---|
307 | cdif khfs est deja w't'_v / heatv(i) = khfs(i) + RETV*zxt*kqfs(i) |
---|
308 | cAM calcule de Ro = paprs(i,1)/Rd zxt |
---|
309 | cAM convention >0 vers le bas ds lmdz |
---|
310 | khfs(i) = - flux_t(i,1)*zxt*Rd / (RCPD*paprs(i,1)) |
---|
311 | kqfs(i) = - flux_q(i,1)*zxt*Rd / (paprs(i,1)) |
---|
312 | cAM verifier que khfs et kqfs sont bien de la forme w'l' |
---|
313 | heatv(i) = khfs(i) + 0.608*zxt*kqfs(i) |
---|
314 | c a comparer aussi aux sorties de clqh : flux_T/RoCp et flux_q/RoLv |
---|
315 | cAM heatv(i) = khfs(i) |
---|
316 | cAM ustar est en entree |
---|
317 | cAM taux = zxu *zxmod*cd_m(i) |
---|
318 | cAM tauy = zxv *zxmod*cd_m(i) |
---|
319 | cAM ustar(i) = SQRT(taux**2+tauy**2) |
---|
320 | cAM ustar(i) = MAX(SQRT(ustar(i)),0.01) |
---|
321 | c Theta et qT du thermique sans exces (interpolin vers surf) |
---|
322 | c chgt de niveau du thermique (jeudi 30/12/1999) |
---|
323 | c (interpolation lineaire avant integration phi_h) |
---|
324 | cAM qT_th(i) = zxqs*beta(i) + 4./z(i,1)*(q(i,1)-zxqs*beta(i)) |
---|
325 | cAM qT_th(i) = max(qT_th(i),q(i,1)) |
---|
326 | qT_th(i) = q2m(i) |
---|
327 | cn The_th restera la Theta du thermique sans exces jusqu'a 2eme calcul |
---|
328 | cn reste a regler convention P) pour Theta |
---|
329 | c The_th(i) = tsol(i) + 4./z(i,1)*(t(i,1)-tsol(i)) |
---|
330 | c - + RLvCp*qT_th(i) |
---|
331 | cAM Th_th(i) = tsol(i) + 4./z(i,1)*(t(i,1)-tsol(i)) |
---|
332 | Th_th(i) = t2m(i) |
---|
333 | ENDDO |
---|
334 | c |
---|
335 | DO i = 1, knon |
---|
336 | rhino(i,1) = 0.0 ! Global Richardson |
---|
337 | check(i) = .TRUE. |
---|
338 | pblh(i) = z(i,1) ! on initialise pblh a l'altitude du 1er niveau |
---|
339 | plcl(i) = 6000. |
---|
340 | c Lambda = -u*^3 / (alpha.g.kvon.<w'Theta'v> |
---|
341 | obklen(i) = -t(i,1)*ustar(i)**3/(RG*vk*heatv(i)) |
---|
342 | trmb1(i) = 0. |
---|
343 | trmb2(i) = 0. |
---|
344 | trmb3(i) = 0. |
---|
345 | ENDDO |
---|
346 | |
---|
347 | C |
---|
348 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
349 | C PBL height calculation: |
---|
350 | C Search for level of pbl. Scan upward until the Richardson number between |
---|
351 | C the first level and the current level exceeds the "critical" value. |
---|
352 | C (bonne idee Nu de separer le Ric et l'exces de temp du thermique) |
---|
353 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
354 | fac = 100.0 |
---|
355 | DO k = 2, isommet |
---|
356 | DO i = 1, knon |
---|
357 | IF (check(i)) THEN |
---|
358 | ! pourquoi / niveau 1 (au lieu du sol) et le terme en u*^2 ? |
---|
359 | ctest zdu2 = (u(i,k)-u(i,1))**2+(v(i,k)-v(i,1))**2+fac*ustar(i)**2 |
---|
360 | zdu2 = u(i,k)**2+v(i,k)**2 |
---|
361 | zdu2 = max(zdu2,1.0e-20) |
---|
362 | c Theta_v environnement |
---|
363 | zthvd=t(i,k)/s(i,k)*(1.+RETV*q(i,k)) |
---|
364 | c |
---|
365 | c therm Theta_v sans exces (avec hypothese fausse de H&B, sinon, |
---|
366 | c passer par Theta_e et virpot) |
---|
367 | c zthvu=t(i,1)/s(i,1)*(1.+RETV*q(i,1)) |
---|
368 | cAM zthvu = Th_th(i)*(1.+RETV*q(i,1)) |
---|
369 | zthvu = Th_th(i)*(1.+RETV*qT_th(i)) |
---|
370 | c Le Ri par Theta_v |
---|
371 | cAM rhino(i,k) = (z(i,k)-z(i,1))*RG*(zthvd-zthvu) |
---|
372 | cAM . /(zdu2*0.5*(zthvd+zthvu)) |
---|
373 | cAM On a nveau de ref a 2m ??? |
---|
374 | rhino(i,k) = (z(i,k)-zref)*RG*(zthvd-zthvu) |
---|
375 | . /(zdu2*0.5*(zthvd+zthvu)) |
---|
376 | c |
---|
377 | IF (rhino(i,k).GE.ricr) THEN |
---|
378 | pblh(i) = z(i,k-1) + (z(i,k-1)-z(i,k)) * |
---|
379 | . (ricr-rhino(i,k-1))/(rhino(i,k-1)-rhino(i,k)) |
---|
380 | c test04 |
---|
381 | pblh(i) = pblh(i) + 100. |
---|
382 | pblT(i) = t(i,k-1) + (t(i,k)-t(i,k-1)) * |
---|
383 | . (pblh(i)-z(i,k-1))/(z(i,k)-z(i,k-1)) |
---|
384 | check(i) = .FALSE. |
---|
385 | ENDIF |
---|
386 | ENDIF |
---|
387 | ENDDO |
---|
388 | ENDDO |
---|
389 | |
---|
390 | C |
---|
391 | C Set pbl height to maximum value where computation exceeds number of |
---|
392 | C layers allowed |
---|
393 | C |
---|
394 | DO i = 1, knon |
---|
395 | if (check(i)) pblh(i) = z(i,isommet) |
---|
396 | ENDDO |
---|
397 | C |
---|
398 | C Improve estimate of pbl height for the unstable points. |
---|
399 | C Find unstable points (sensible heat flux is upward): |
---|
400 | C |
---|
401 | DO i = 1, knon |
---|
402 | IF (heatv(i) .GT. 0.) THEN |
---|
403 | unstbl(i) = .TRUE. |
---|
404 | check(i) = .TRUE. |
---|
405 | ELSE |
---|
406 | unstbl(i) = .FALSE. |
---|
407 | check(i) = .FALSE. |
---|
408 | ENDIF |
---|
409 | ENDDO |
---|
410 | C |
---|
411 | C For the unstable case, compute velocity scale and the |
---|
412 | C convective temperature excess: |
---|
413 | C |
---|
414 | DO i = 1, knon |
---|
415 | IF (check(i)) THEN |
---|
416 | phiminv(i) = (1.-binm*pblh(i)/obklen(i))**onet |
---|
417 | c *************************************************** |
---|
418 | c Wm ? et W* ? c'est la formule pour z/h < .1 |
---|
419 | c ;Calcul de w* ;; |
---|
420 | c ;;;;;;;;;;;;;;;; |
---|
421 | c w_star=((g/tcls)*fcsv*z(ind))^(1/3.) [ou prendre la premiere approx de h) |
---|
422 | c ;; CALCUL DE wm ;; |
---|
423 | c ;;;;;;;;;;;;;;;;;; |
---|
424 | c ; Ici on considerera que l'on est dans la couche de surf jusqu'a 100m |
---|
425 | c ; On prend svt couche de surface=0.1*h mais on ne connait pas h |
---|
426 | c ;;;;;;;;;;;Dans la couche de surface |
---|
427 | c if (z(ind) le 20) then begin |
---|
428 | c Phim=(1.-15.*(z(ind)/L))^(-1/3.) |
---|
429 | c wm=u_star/Phim |
---|
430 | c ;;;;;;;;;;;En dehors de la couche de surface |
---|
431 | c endif else if (z(ind) gt 20) then begin |
---|
432 | c wm=(u_star^3+c1*w_star^3)^(1/3.) |
---|
433 | c endif |
---|
434 | c *************************************************** |
---|
435 | wm(i)= ustar(i)*phiminv(i) |
---|
436 | c====================================================================== |
---|
437 | cvaleurs de Dominique Lambert de la campagne SEMAPHORE : |
---|
438 | c <T'^2> = 100.T*^2; <q'^2> = 20.q*^2 a 10m |
---|
439 | c <Tv'^2> = (1+1.2q).100.T* + 1.2Tv.sqrt(20*100).T*.q* + (.608*Tv)^2*20.q*^2; |
---|
440 | c et dTetavS = sqrt(<Tv'^2>) ainsi calculee. |
---|
441 | c avec : T*=<w'T'>_s/w* et q*=<w'q'>/w* |
---|
442 | c !!! on peut donc utiliser w* pour les fluctuations <-> Lambert |
---|
443 | c(leur corellation pourrait dependre de beta par ex) |
---|
444 | c if fcsv(i,j) gt 0 then begin |
---|
445 | c dTetavs=b1*(1.+2.*.608*q_10(i,j))*(fcs(i,j)/wm(i,j))^2+$ |
---|
446 | c (.608*Thetav_10(i,j))^2*b2*(xle(i,j)/wm(i,j))^2+$ |
---|
447 | c 2.*.608*thetav_10(i,j)*sqrt(b1*b2)*(xle(i,j)/wm(i,j))*(fcs(i,j)/wm(i,j)) |
---|
448 | c dqs=b2*(xle(i,j)/wm(i,j))^2 |
---|
449 | c theta_s(i,j)=thetav_10(i,j)+sqrt(dTetavs) |
---|
450 | c q_s(i,j)=q_10(i,j)+sqrt(dqs) |
---|
451 | c endif else begin |
---|
452 | c Theta_s(i,j)=thetav_10(i,j) |
---|
453 | c q_s(i,j)=q_10(i,j) |
---|
454 | c endelse |
---|
455 | c====================================================================== |
---|
456 | c |
---|
457 | cHBTM therm(i) = heatv(i)*fak/wm(i) |
---|
458 | c forme Mathieu : |
---|
459 | q_star = kqfs(i)/wm(i) |
---|
460 | t_star = khfs(i)/wm(i) |
---|
461 | cIM 091204 BEG |
---|
462 | IF(1.EQ.0) THEN |
---|
463 | IF(t_star.LT.0..OR.q_star.LT.0.) THEN |
---|
464 | print*,'i t_star q_star khfs kqfs wm',i,t_star,q_star, |
---|
465 | $ khfs(i),kqfs(i),wm(i) |
---|
466 | ENDIF |
---|
467 | ENDIF |
---|
468 | cIM 091204 END |
---|
469 | cAM Nveau cde ref 2m => |
---|
470 | cAM therm(i) = sqrt( b1*(1.+2.*RETV*q(i,1))*t_star**2 |
---|
471 | cAM + + (RETV*T(i,1))**2*b2*q_star**2 |
---|
472 | cAM + + 2.*RETV*T(i,1)*b212*q_star*t_star |
---|
473 | cAM + ) |
---|
474 | cIM 091204 BEG |
---|
475 | a1=b1*(1.+2.*RETV*qT_th(i))*t_star**2 |
---|
476 | a2=(RETV*Th_th(i))**2*b2*q_star**2 |
---|
477 | a3=2.*RETV*Th_th(i)*b212*q_star*t_star |
---|
478 | aa=a1+a2+a3 |
---|
479 | IF(1.EQ.0) THEN |
---|
480 | IF (aa.LT.0.) THEN |
---|
481 | print*,'i a1 a2 a3 aa',i,a1,a2,a3,aa |
---|
482 | print*,'i qT_th Th_th t_star q_star RETV b1 b2 b212', |
---|
483 | $ i,qT_th(i),Th_th(i),t_star,q_star,RETV,b1,b2,b212 |
---|
484 | ENDIF |
---|
485 | ENDIF |
---|
486 | cIM 091204 END |
---|
487 | therm(i) = sqrt( b1*(1.+2.*RETV*qT_th(i))*t_star**2 |
---|
488 | + + (RETV*Th_th(i))**2*b2*q_star**2 |
---|
489 | cIM 101204 + + 2.*RETV*Th_th(i)*b212*q_star*t_star |
---|
490 | + + max(0.,2.*RETV*Th_th(i)*b212*q_star*t_star) |
---|
491 | + ) |
---|
492 | c |
---|
493 | c Theta et qT du thermique (forme H&B) avec exces |
---|
494 | c (attention, on ajoute therm(i) qui est virtuelle ...) |
---|
495 | c pourquoi pas sqrt(b1)*t_star ? |
---|
496 | c dqs = b2sr*kqfs(i)/wm(i) |
---|
497 | qT_th(i) = qT_th(i) + b2sr*q_star |
---|
498 | cnew on differre le calcul de Theta_e |
---|
499 | c The_th(i) = The_th(i) + therm(i) + RLvCp*qT_th(i) |
---|
500 | c ou: The_th(i) = The_th(i) + sqrt(b1)*khfs(i)/wm(i) + RLvCp*qT_th(i) |
---|
501 | rhino(i,1) = 0.0 |
---|
502 | ENDIF |
---|
503 | ENDDO |
---|
504 | C |
---|
505 | c +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
506 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
507 | C ++ Improve pblh estimate for unstable conditions using the +++++++ |
---|
508 | C ++ convective temperature excess : +++++++ |
---|
509 | c ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
510 | c +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
---|
511 | C |
---|
512 | DO k = 2, isommet |
---|
513 | DO i = 1, knon |
---|
514 | IF (check(i)) THEN |
---|
515 | ctest zdu2 = (u(i,k)-u(i,1))**2+(v(i,k)-v(i,1))**2+fac*ustar(i)**2 |
---|
516 | zdu2 = u(i,k)**2+v(i,k)**2 |
---|
517 | zdu2 = max(zdu2,1.0e-20) |
---|
518 | c Theta_v environnement |
---|
519 | zthvd=t(i,k)/s(i,k)*(1.+RETV*q(i,k)) |
---|
520 | c |
---|
521 | c et therm Theta_v (avec hypothese de constance de H&B, |
---|
522 | c zthvu=(t(i,1)+therm(i))/s(i,1)*(1.+RETV*q(i,1)) |
---|
523 | zthvu = Th_th(i)*(1.+RETV*qT_th(i)) + therm(i) |
---|
524 | |
---|
525 | c |
---|
526 | c Le Ri par Theta_v |
---|
527 | cAM Niveau de ref 2m |
---|
528 | cAM rhino(i,k) = (z(i,k)-z(i,1))*RG*(zthvd-zthvu) |
---|
529 | cAM . /(zdu2*0.5*(zthvd+zthvu)) |
---|
530 | rhino(i,k) = (z(i,k)-zref)*RG*(zthvd-zthvu) |
---|
531 | . /(zdu2*0.5*(zthvd+zthvu)) |
---|
532 | c |
---|
533 | c |
---|
534 | IF (rhino(i,k).GE.ricr) THEN |
---|
535 | pblh(i) = z(i,k-1) + (z(i,k-1)-z(i,k)) * |
---|
536 | . (ricr-rhino(i,k-1))/(rhino(i,k-1)-rhino(i,k)) |
---|
537 | c test04 |
---|
538 | pblh(i) = pblh(i) + 100. |
---|
539 | pblT(i) = t(i,k-1) + (t(i,k)-t(i,k-1)) * |
---|
540 | . (pblh(i)-z(i,k-1))/(z(i,k)-z(i,k-1)) |
---|
541 | check(i) = .FALSE. |
---|
542 | cIM 170305 BEG |
---|
543 | IF(1.EQ.0) THEN |
---|
544 | c debug print -120;34 -34- 58 et 0;26 wamp |
---|
545 | if (i.eq.950.or.i.eq.192.or.i.eq.624.or.i.eq.118) then |
---|
546 | print*,' i,Th_th,Therm,qT :',i,Th_th(i),therm(i),qT_th(i) |
---|
547 | q_star = kqfs(i)/wm(i) |
---|
548 | t_star = khfs(i)/wm(i) |
---|
549 | print*,'q* t*, b1,b2,b212 ',q_star,t_star |
---|
550 | - , b1*(1.+2.*RETV*qT_th(i))*t_star**2 |
---|
551 | - , (RETV*Th_th(i))**2*b2*q_star**2 |
---|
552 | - , 2.*RETV*Th_th(i)*b212*q_star*t_star |
---|
553 | print*,'zdu2 ,100.*ustar(i)**2',zdu2 ,fac*ustar(i)**2 |
---|
554 | endif |
---|
555 | ENDIF !(1.EQ.0) THEN |
---|
556 | cIM 170305 END |
---|
557 | c q_star = kqfs(i)/wm(i) |
---|
558 | c t_star = khfs(i)/wm(i) |
---|
559 | c trmb1(i) = b1*(1.+2.*RETV*q(i,1))*t_star**2 |
---|
560 | c trmb2(i) = (RETV*T(i,1))**2*b2*q_star**2 |
---|
561 | c Omega now trmb3(i) = 2.*RETV*T(i,1)*b212*q_star*t_star |
---|
562 | ENDIF |
---|
563 | ENDIF |
---|
564 | ENDDO |
---|
565 | ENDDO |
---|
566 | C |
---|
567 | C Set pbl height to maximum value where computation exceeds number of |
---|
568 | C layers allowed |
---|
569 | C |
---|
570 | DO i = 1, knon |
---|
571 | if (check(i)) pblh(i) = z(i,isommet) |
---|
572 | ENDDO |
---|
573 | C |
---|
574 | C PBL height must be greater than some minimum mechanical mixing depth |
---|
575 | C Several investigators have proposed minimum mechanical mixing depth |
---|
576 | C relationships as a function of the local friction velocity, u*. We |
---|
577 | C make use of a linear relationship of the form h = c u* where c=700. |
---|
578 | C The scaling arguments that give rise to this relationship most often |
---|
579 | C represent the coefficient c as some constant over the local coriolis |
---|
580 | C parameter. Here we make use of the experimental results of Koracin |
---|
581 | C and Berkowicz (1988) [BLM, Vol 43] for wich they recommend 0.07/f |
---|
582 | C where f was evaluated at 39.5 N and 52 N. Thus we use a typical mid |
---|
583 | C latitude value for f so that c = 0.07/f = 700. |
---|
584 | C |
---|
585 | DO i = 1, knon |
---|
586 | pblmin = 700.0*ustar(i) |
---|
587 | pblh(i) = MAX(pblh(i),pblmin) |
---|
588 | c par exemple : |
---|
589 | pblT(i) = t(i,2) + (t(i,3)-t(i,2)) * |
---|
590 | . (pblh(i)-z(i,2))/(z(i,3)-z(i,2)) |
---|
591 | ENDDO |
---|
592 | |
---|
593 | C ******************************************************************** |
---|
594 | C pblh is now available; do preparation for diffusivity calculation : |
---|
595 | C ******************************************************************** |
---|
596 | DO i = 1, knon |
---|
597 | check(i) = .TRUE. |
---|
598 | Zsat(i) = .FALSE. |
---|
599 | c omegafl utilise pour prolongement CAPE |
---|
600 | omegafl(i) = .FALSE. |
---|
601 | Cape(i) = 0. |
---|
602 | Kape(i) = 0. |
---|
603 | EauLiq(i) = 0. |
---|
604 | CTEI(i) = 0. |
---|
605 | pblk(i) = 0.0 |
---|
606 | fak1(i) = ustar(i)*pblh(i)*vk |
---|
607 | C |
---|
608 | C Do additional preparation for unstable cases only, set temperature |
---|
609 | C and moisture perturbations depending on stability. |
---|
610 | C *** Rq: les formule sont prises dans leur forme CS *** |
---|
611 | IF (unstbl(i)) THEN |
---|
612 | cAM Niveau de ref du thermique |
---|
613 | cAM zxt=(t(i,1)-z(i,1)*0.5*RG/RCPD/(1.+RVTMP2*q(i,1))) |
---|
614 | cAM . *(1.+RETV*q(i,1)) |
---|
615 | zxt=(Th_th(i)-zref*0.5*RG/RCPD/(1.+RVTMP2*qT_th(i))) |
---|
616 | . *(1.+RETV*qT_th(i)) |
---|
617 | phiminv(i) = (1. - binm*pblh(i)/obklen(i))**onet |
---|
618 | phihinv(i) = sqrt(1. - binh*pblh(i)/obklen(i)) |
---|
619 | wm(i) = ustar(i)*phiminv(i) |
---|
620 | fak2(i) = wm(i)*pblh(i)*vk |
---|
621 | wstr(i) = (heatv(i)*RG*pblh(i)/zxt)**onet |
---|
622 | fak3(i) = fakn*wstr(i)/wm(i) |
---|
623 | ENDIF |
---|
624 | c Computes Theta_e for thermal (all cases : to be modified) |
---|
625 | c attention ajout therm(i) = virtuelle |
---|
626 | The_th(i) = Th_th(i) + therm(i) + RLvCp*qT_th(i) |
---|
627 | c ou: The_th(i) = Th_th(i) + sqrt(b1)*khfs(i)/wm(i) + RLvCp*qT_th(i) |
---|
628 | ENDDO |
---|
629 | |
---|
630 | C Main level loop to compute the diffusivities and |
---|
631 | C counter-gradient terms: |
---|
632 | C |
---|
633 | DO 1000 k = 2, isommet |
---|
634 | C |
---|
635 | C Find levels within boundary layer: |
---|
636 | C |
---|
637 | DO i = 1, knon |
---|
638 | unslev(i) = .FALSE. |
---|
639 | stblev(i) = .FALSE. |
---|
640 | zm(i) = z(i,k-1) |
---|
641 | zp(i) = z(i,k) |
---|
642 | IF (zkmin.EQ.0.0 .AND. zp(i).GT.pblh(i)) zp(i) = pblh(i) |
---|
643 | IF (zm(i) .LT. pblh(i)) THEN |
---|
644 | zmzp = 0.5*(zm(i) + zp(i)) |
---|
645 | C debug |
---|
646 | c if (i.EQ.1864) then |
---|
647 | c print*,'i,pblh(1864),obklen(1864)',i,pblh(i),obklen(i) |
---|
648 | c endif |
---|
649 | |
---|
650 | zh(i) = zmzp/pblh(i) |
---|
651 | zl(i) = zmzp/obklen(i) |
---|
652 | zzh(i) = 0. |
---|
653 | IF (zh(i).LE.1.0) zzh(i) = (1. - zh(i))**2 |
---|
654 | C |
---|
655 | C stblev for points zm < plbh and stable and neutral |
---|
656 | C unslev for points zm < plbh and unstable |
---|
657 | C |
---|
658 | IF (unstbl(i)) THEN |
---|
659 | unslev(i) = .TRUE. |
---|
660 | ELSE |
---|
661 | stblev(i) = .TRUE. |
---|
662 | ENDIF |
---|
663 | ENDIF |
---|
664 | ENDDO |
---|
665 | c print*,'fin calcul niveaux' |
---|
666 | C |
---|
667 | C Stable and neutral points; set diffusivities; counter-gradient |
---|
668 | C terms zero for stable case: |
---|
669 | C |
---|
670 | DO i = 1, knon |
---|
671 | IF (stblev(i)) THEN |
---|
672 | IF (zl(i).LE.1.) THEN |
---|
673 | pblk(i) = fak1(i)*zh(i)*zzh(i)/(1. + betas*zl(i)) |
---|
674 | ELSE |
---|
675 | pblk(i) = fak1(i)*zh(i)*zzh(i)/(betas + zl(i)) |
---|
676 | ENDIF |
---|
677 | c pcfm(i,k) = pblk(i) |
---|
678 | c pcfh(i,k) = pcfm(i,k) |
---|
679 | ENDIF |
---|
680 | ENDDO |
---|
681 | C |
---|
682 | C unssrf, unstable within surface layer of pbl |
---|
683 | C unsout, unstable within outer layer of pbl |
---|
684 | C |
---|
685 | DO i = 1, knon |
---|
686 | unssrf(i) = .FALSE. |
---|
687 | unsout(i) = .FALSE. |
---|
688 | IF (unslev(i)) THEN |
---|
689 | IF (zh(i).lt.sffrac) THEN |
---|
690 | unssrf(i) = .TRUE. |
---|
691 | ELSE |
---|
692 | unsout(i) = .TRUE. |
---|
693 | ENDIF |
---|
694 | ENDIF |
---|
695 | ENDDO |
---|
696 | C |
---|
697 | C Unstable for surface layer; counter-gradient terms zero |
---|
698 | C |
---|
699 | DO i = 1, knon |
---|
700 | IF (unssrf(i)) THEN |
---|
701 | term = (1. - betam*zl(i))**onet |
---|
702 | pblk(i) = fak1(i)*zh(i)*zzh(i)*term |
---|
703 | pr(i) = term/sqrt(1. - betah*zl(i)) |
---|
704 | ENDIF |
---|
705 | ENDDO |
---|
706 | c print*,'fin counter-gradient terms zero' |
---|
707 | C |
---|
708 | C Unstable for outer layer; counter-gradient terms non-zero: |
---|
709 | C |
---|
710 | DO i = 1, knon |
---|
711 | IF (unsout(i)) THEN |
---|
712 | pblk(i) = fak2(i)*zh(i)*zzh(i) |
---|
713 | c cgs(i,k) = fak3(i)/(pblh(i)*wm(i)) |
---|
714 | c cgh(i,k) = khfs(i)*cgs(i,k) |
---|
715 | pr(i) = phiminv(i)/phihinv(i) + ccon*fak3(i)/fak |
---|
716 | c cgq(i,k) = kqfs(i)*cgs(i,k) |
---|
717 | ENDIF |
---|
718 | ENDDO |
---|
719 | c print*,'fin counter-gradient terms non zero' |
---|
720 | C |
---|
721 | C For all unstable layers, compute diffusivities and ctrgrad ter m |
---|
722 | C |
---|
723 | c DO i = 1, knon |
---|
724 | c IF (unslev(i)) THEN |
---|
725 | c pcfm(i,k) = pblk(i) |
---|
726 | c pcfh(i,k) = pblk(i)/pr(i) |
---|
727 | c etc cf original |
---|
728 | c ENDIF |
---|
729 | c ENDDO |
---|
730 | C |
---|
731 | C For all layers, compute integral info and CTEI |
---|
732 | C |
---|
733 | DO i = 1, knon |
---|
734 | if (check(i).or.omegafl(i)) then |
---|
735 | if (.not.Zsat(i)) then |
---|
736 | c Th2 = The_th(i) - RLvCp*qT_th(i) |
---|
737 | Th2 = Th_th(i) |
---|
738 | T2 = Th2*s(i,k) |
---|
739 | c thermodyn functions |
---|
740 | zdelta=MAX(0.,SIGN(1.,RTT-T2)) |
---|
741 | qqsat= r2es * FOEEW(T2,zdelta)/pplay(i,k) |
---|
742 | qqsat=MIN(0.5,qqsat) |
---|
743 | zcor=1./(1.-retv*qqsat) |
---|
744 | qqsat=qqsat*zcor |
---|
745 | c |
---|
746 | if (qqsat.lt.qT_th(i)) then |
---|
747 | c on calcule lcl |
---|
748 | if (k.eq.2) then |
---|
749 | plcl(i) = z(i,k) |
---|
750 | else |
---|
751 | plcl(i) = z(i,k-1) + (z(i,k-1)-z(i,k)) * |
---|
752 | . (qT_th(i)-qsatbef(i))/(qsatbef(i)-qqsat) |
---|
753 | endif |
---|
754 | Zsat(i) = .true. |
---|
755 | Tbef(i) = T2 |
---|
756 | endif |
---|
757 | c |
---|
758 | endif |
---|
759 | qsatbef(i) = qqsat |
---|
760 | camn ???? cette ligne a deja ete faite normalement ? |
---|
761 | endif |
---|
762 | c print*,'hbtm2 i,k=',i,k |
---|
763 | ENDDO |
---|
764 | 1000 continue ! end of level loop |
---|
765 | cIM 170305 BEG |
---|
766 | IF(1.EQ.0) THEN |
---|
767 | print*,'hbtm2 ok' |
---|
768 | ENDIF !(1.EQ.0) THEN |
---|
769 | cIM 170305 END |
---|
770 | RETURN |
---|
771 | END |
---|