1 | MODULE pbl_parameters_mod |
---|
2 | |
---|
3 | !======================================================================================================================! |
---|
4 | ! Subject: |
---|
5 | !--------- |
---|
6 | ! Module used to compute the friction velocity, temperature, and monin_obukhov length from temperature, wind field and thermals outputs. |
---|
7 | ! Also interpolate theta and u;v in the surface layer based on the Monin Obukhov theory |
---|
8 | ! These are diagnostics only and do not influence the code. |
---|
9 | !----------------------------------------------------------------------------------------------------------------------! |
---|
10 | ! Reference: |
---|
11 | !----------- |
---|
12 | ! Colaïtis, A., Spiga, A., Hourdin, F., Rio, C., Forget, F., and Millour, E. (2013), A thermal plume model for the Martian convective boundary layer, J. Geophys. Res. Planets, 118, 1468–1487, doi:10.1002/jgre.20104. |
---|
13 | ! |
---|
14 | !======================================================================================================================! |
---|
15 | |
---|
16 | |
---|
17 | IMPLICIT NONE |
---|
18 | |
---|
19 | CONTAINS |
---|
20 | |
---|
21 | SUBROUTINE pbl_parameters(ngrid,nlay,ps,pplay,pz0,pg,zzlay,zzlev,pu,pv,wstar_in,hfmax,zmax,tke,pts,ph,pqvap,pqsurf,mumean,z_out,n_out, & |
---|
22 | T_out,u_out,ustar,tstar,vhf,vvv) |
---|
23 | |
---|
24 | USE comcstfi_h |
---|
25 | use write_output_mod, only: write_output |
---|
26 | use turb_mod, only: turb_resolved |
---|
27 | use lmdz_atke_turbulence_ini, only : smmin, ric, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi |
---|
28 | use watersat_mod, only: watersat |
---|
29 | use paleoclimate_mod, only: include_waterbuoyancy |
---|
30 | |
---|
31 | IMPLICIT NONE |
---|
32 | |
---|
33 | !======================================================================= |
---|
34 | ! |
---|
35 | ! Anlysis of the PBL from input temperature, wind field and thermals outputs: compute the friction velocity, temperature, and monin_obukhov length |
---|
36 | ! and interpolate the potential temperature and winds in the surface layer using the Monin Obukhov theory |
---|
37 | ! |
---|
38 | ! ------- |
---|
39 | ! |
---|
40 | ! Author: Arnaud Colaitis 09/01/12; adapted in F90 by Lucas Lange 12/2023 |
---|
41 | ! ------- |
---|
42 | ! |
---|
43 | ! Arguments: |
---|
44 | ! ---------- |
---|
45 | ! |
---|
46 | ! inputs: |
---|
47 | ! ------ |
---|
48 | ! ngrid size of the horizontal grid |
---|
49 | ! nlay size of the vertical grid |
---|
50 | ! ps(ngrid) surface pressure (Pa) |
---|
51 | ! pplay(ngrid,nlay) pressure levels |
---|
52 | ! pz0(ngrid) surface roughness length (m) |
---|
53 | ! pg gravity (m s -2) |
---|
54 | ! zzlay(ngrid,nlay) height of mid-layers (m) |
---|
55 | ! zzlev(ngrid,nlay+1) height of layers interfaces (m) |
---|
56 | ! pu(ngrid,nlay) u component of the wind (m/s) |
---|
57 | ! pv(ngrid,nlay) v component of the wind (m/s) |
---|
58 | ! wstar_in(ngrid) free convection velocity in PBL (m/s) |
---|
59 | ! hfmax(ngrid) maximum vertical turbulent heat flux in thermals (W/m^2) |
---|
60 | ! zmax(ngrid) height reached by the thermals (pbl height) (m) |
---|
61 | ! tke(ngrid,nlay+1) turbulent kinetic energy (J/kg) |
---|
62 | ! pts(ngrid) surface temperature (K) |
---|
63 | ! ph(ngrid,nlay) potential temperature T*(p/ps)^kappa (k) |
---|
64 | ! z_out(n_out) heights of interpolation (m) |
---|
65 | ! n_out number of points for interpolation (1) |
---|
66 | ! |
---|
67 | ! outputs: |
---|
68 | ! ------ |
---|
69 | ! |
---|
70 | ! T_out(ngrid,n_out) interpolated potential temperature on z_out (K) |
---|
71 | ! u_out(ngrid,n_out) interpolated winds on z_out (m/s) |
---|
72 | ! ustar(ngrid) friction velocity (m/s) |
---|
73 | ! tstar(ngrid) friction temperature (K) |
---|
74 | ! |
---|
75 | ! |
---|
76 | !======================================================================= |
---|
77 | ! |
---|
78 | !----------------------------------------------------------------------- |
---|
79 | ! Declarations: |
---|
80 | ! ------------- |
---|
81 | |
---|
82 | #include "callkeys.h" |
---|
83 | |
---|
84 | ! Arguments: |
---|
85 | ! ---------- |
---|
86 | |
---|
87 | INTEGER, INTENT(IN) :: ngrid,nlay,n_out ! size of the horizontal and vertical grid, interpolated altitudes for the surface layer |
---|
88 | REAL, INTENT(IN) :: pz0(ngrid) ! surface roughness |
---|
89 | REAL, INTENT(IN) :: ps(ngrid),pplay(ngrid,nlay) ! surface pressure and pressure levels |
---|
90 | REAL, INTENT(IN) :: pg,zzlay(ngrid,nlay),zzlev(ngrid,nlay) ! surface gravity, altitude of the interface and mid-layers |
---|
91 | REAL, INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay) ! winds |
---|
92 | REAL, INTENT(IN) :: wstar_in(ngrid),hfmax(ngrid),zmax(ngrid) ! free convection velocity , vertical turbulent heat, pbl height |
---|
93 | REAL, INTENT(IN) :: tke(ngrid,nlay+1) ! Turbulent kinetic energy |
---|
94 | REAL, INTENT(IN) :: pts(ngrid),ph(ngrid,nlay) ! surface temperature, potentiel temperature |
---|
95 | REAL, INTENT(IN) :: z_out(n_out) ! altitudes of the interpolation in the surface layer |
---|
96 | REAL, INTENT(IN) :: mumean(ngrid) ! Molecular mass of the atmosphere (kg/mol) |
---|
97 | REAL, INTENT(IN) :: pqvap(ngrid,nlay) ! Water vapor mixing ratio (kg/kg) to account for vapor flottability |
---|
98 | REAL, INTENT(IN) :: pqsurf(ngrid) ! Water ice frost on the surface (kg/m^2) to account for vapor flottability |
---|
99 | |
---|
100 | ! Outputs: |
---|
101 | ! -------- |
---|
102 | |
---|
103 | REAL, INTENT(OUT) :: T_out(ngrid,n_out),u_out(ngrid,n_out) ! Temperature and wind of the interpolation in the surface layer |
---|
104 | REAL, INTENT(OUT) :: ustar(ngrid), tstar(ngrid) ! Friction velocity and temperature |
---|
105 | |
---|
106 | ! Local: |
---|
107 | ! ------ |
---|
108 | |
---|
109 | INTEGER ig,k,n ! Loop variables |
---|
110 | INTEGER ii(1) ! Index to compute the pbl mixed layer temperature |
---|
111 | REAL karman,nu ! Von Karman constant and fluid kinematic viscosity |
---|
112 | DATA karman,nu/.41,0.001/ |
---|
113 | |
---|
114 | !$OMP THREADPRIVATE(karman,nu) |
---|
115 | |
---|
116 | SAVE karman,nu |
---|
117 | |
---|
118 | REAL zout ! altitude to interpolate (local variable during the loop on z_out) (m) |
---|
119 | REAL rib(ngrid) ! Bulk Richardson number (1) |
---|
120 | REAL fm(ngrid) ! stability function for momentum (1) |
---|
121 | REAL fh(ngrid) ! stability function for heat (1) |
---|
122 | REAL z1z0,z1z0t ! ratios z1/z0 and z1/z0T (1) |
---|
123 | ! Formalism for stability functions fm= 1/phim^2; fh = 1/(phim*phih) where |
---|
124 | ! phim = 1+betam*zeta or (1-bm*zeta)**am |
---|
125 | ! phih = alphah + betah*zeta or alphah(1.-bh*zeta)**ah |
---|
126 | ! ah and am are assumed to be -0.25 and -0.5 respectively |
---|
127 | ! lambda is an intermediate variable to simplify the expression |
---|
128 | REAL betam, betah, alphah, bm, bh, lambda |
---|
129 | |
---|
130 | REAL cdn(ngrid),chn(ngrid) ! neutral momentum and heat drag coefficient (1) |
---|
131 | REAL pz0t ! initial thermal roughness length z0t for the iterative scheme (m) |
---|
132 | REAL ric_colaitis ! critical richardson number (1) |
---|
133 | REAL reynolds(ngrid) ! Reynolds number (1) |
---|
134 | REAL prandtl(ngrid) ! Prandtl number (1) |
---|
135 | REAL pz0tcomp(ngrid) ! Computed z0t during the iterative scheme (m) |
---|
136 | REAL ite ! Number of iteration (1) |
---|
137 | REAL residual ! Residual between pz0t and pz0tcomp (m) |
---|
138 | REAL itemax ! Maximum number of iteration (1) |
---|
139 | REAL tol_iter ! Tolerance for the residual to ensure convergence (1) |
---|
140 | REAL pcdv(ngrid),pcdh(ngrid) ! momentum and heat drag coefficient (1) |
---|
141 | REAL zu2(ngrid) ! Near surface winds (m/s) |
---|
142 | REAL x(ngrid) ! z/zi (1) |
---|
143 | REAL dvhf(ngrid),dvvv(ngrid) ! dimensionless vertical heat flux and dimensionless vertical velocity variance |
---|
144 | REAL vhf(ngrid),vvv(ngrid) ! vertical heat flux (W/m^2) and vertical velocity variance (m/s) |
---|
145 | REAL Teta_out(ngrid,n_out) ! Temporary variable to compute interpolated potential temperature (K) |
---|
146 | REAL sm ! Stability function computed with the ATKE scheme |
---|
147 | REAL f_ri_cd_min ! minimum of the stability function with the ATKE scheme |
---|
148 | REAL ric_4interp ! critical richardson number which is either the one from Colaitis et al. (2013) or ATKE (1) |
---|
149 | |
---|
150 | |
---|
151 | REAL tsurf_v(ngrid) ! Virtual surface temperature, accounting for vapor flottability |
---|
152 | REAL temp_v(ngrid) ! Potential virtual air temperature in the frist layer, accounting for vapor flottability |
---|
153 | REAL mu_h2o ! Molecular mass of water vapor (kg/mol) |
---|
154 | REAL tol_frost ! Tolerance to consider the effect of frost on the surface |
---|
155 | REAL qsat(ngrid) ! saturated mixing ratio of water vapor |
---|
156 | |
---|
157 | |
---|
158 | ! Code: |
---|
159 | ! -------- |
---|
160 | |
---|
161 | !------------------------------------------------------------------------ |
---|
162 | !------------------------------------------------------------------------ |
---|
163 | ! PART I : RICHARDSON/REYNOLDS/THERMAL_ROUGHNESS/STABILITY_COEFFICIENTS |
---|
164 | !------------------------------------------------------------------------ |
---|
165 | !------------------------------------------------------------------------ |
---|
166 | |
---|
167 | |
---|
168 | |
---|
169 | ! Initialisation : |
---|
170 | |
---|
171 | ustar(:)=0. |
---|
172 | tstar(:)=0. |
---|
173 | reynolds(:)=0. |
---|
174 | pz0t = 0. |
---|
175 | pz0tcomp(:) = 0.1*pz0(:) |
---|
176 | rib(:)=0. |
---|
177 | pcdv(:)=0. |
---|
178 | pcdh(:)=0. |
---|
179 | itemax= 10 |
---|
180 | tol_iter = 0.01 |
---|
181 | f_ri_cd_min = 0.01 |
---|
182 | mu_h2o = 18e-3 |
---|
183 | tol_frost = 1e-4 |
---|
184 | tsurf_v(:) = 0. |
---|
185 | temp_v(:) = 0. |
---|
186 | |
---|
187 | ! this formulation assumes alphah=1., implying betah=betam |
---|
188 | ! We use Dyer et al. parameters, as they cover a broad range of Richardson numbers : |
---|
189 | bm = 16. |
---|
190 | bh = 16. |
---|
191 | alphah = 1. |
---|
192 | betam = 5. |
---|
193 | betah = 5. |
---|
194 | lambda = (sqrt(bh/bm))/alphah |
---|
195 | ric_colaitis = betah/(betam**2) |
---|
196 | |
---|
197 | if(callatke) then |
---|
198 | ric_4interp = ric |
---|
199 | else |
---|
200 | ric_4interp = ric_colaitis |
---|
201 | endif |
---|
202 | |
---|
203 | |
---|
204 | |
---|
205 | IF(include_waterbuoyancy) then |
---|
206 | temp_v(:) = ph(:,1)*(1.+pqvap(:,1)/(mu_h2o/mumean(:)))/(1.+pqvap(:,1)) |
---|
207 | DO ig = 1,ngrid |
---|
208 | IF(pqsurf(ig).gt.tol_frost) then |
---|
209 | call watersat(1,pts(ig),pplay(ig,1),qsat(ig)) |
---|
210 | tsurf_v(ig) = pts(ig)*(1.+qsat(ig)/(mu_h2o/mumean(ig)))/(1.+qsat(ig)) |
---|
211 | ELSE |
---|
212 | tsurf_v(ig) = pts(ig)*(1.+pqvap(ig,1)/(mu_h2o/mumean(ig)))/(1.+pqvap(ig,1)) |
---|
213 | ENDIF |
---|
214 | ENDDO |
---|
215 | ELSE |
---|
216 | tsurf_v(:) = pts(:) |
---|
217 | temp_v(:) = ph(:,1) |
---|
218 | ENDIF |
---|
219 | |
---|
220 | DO ig=1,ngrid |
---|
221 | ite = 0. |
---|
222 | residual = abs(pz0tcomp(ig)-pz0t) |
---|
223 | z1z0 = zzlay(ig,1)/pz0(ig) |
---|
224 | cdn(ig) = karman/log(z1z0) |
---|
225 | cdn(ig) = cdn(ig)*cdn(ig) |
---|
226 | zu2(ig) = pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1) + (log(1.+0.7*wstar_in(ig) + 2.3*wstar_in(ig)**2))**2 ! near surface winds, accounting for buyoncy |
---|
227 | ! IF(turb_resolved) zu2(ig)=MAX(zu2(ig),1.) |
---|
228 | DO WHILE((residual .gt. tol_iter*pz0(ig)) .and. (ite .lt. itemax)) |
---|
229 | ! Computations of z0T; iterated until z0T converges |
---|
230 | pz0t = pz0tcomp(ig) |
---|
231 | z1z0t=zzlay(ig,1)/pz0t |
---|
232 | chn(ig) = cdn(ig)*log(z1z0)/log(z1z0t) |
---|
233 | IF ((pu(ig,1) .ne. 0.) .or. (pv(ig,1) .ne. 0.)) then |
---|
234 | ! Richardson number formulation proposed by D.E. England et al. (1995) |
---|
235 | rib(ig) = (pg/tsurf_v(ig))*sqrt(zzlay(ig,1)*pz0(ig))*(((log(zzlay(ig,1)/pz0(ig)))**2)/(log(zzlay(ig,1)/pz0t)))*(temp_v(ig)-tsurf_v(ig))/zu2(ig) |
---|
236 | ELSE |
---|
237 | print*,'warning, infinite Richardson at surface' |
---|
238 | print*,pu(ig,1),pv(ig,1) |
---|
239 | rib(ig) = ric_colaitis |
---|
240 | ENDIF |
---|
241 | |
---|
242 | ! Compute the stability functions fm; fh depending on the stability of the surface layer |
---|
243 | |
---|
244 | IF(callatke) THEN |
---|
245 | ! Computation following parameterizaiton from ATKE |
---|
246 | IF(rib(ig) .gt. 0) THEN |
---|
247 | ! STABLE BOUNDARY LAYER : |
---|
248 | sm = MAX(smmin,cn*(1.-rib(ig)/ric)) |
---|
249 | ! prandlt expression from venayagamoorthy and stretch 2010, Li et al 2019 |
---|
250 | prandtl(ig) = pr_neut*exp(-pr_slope/pr_neut*rib(ig)+rib(ig)/pr_neut) + rib(ig) * pr_slope |
---|
251 | fm(ig) = max((sm**(3./2.) * sqrt(cepsilon) * (1 - rib(ig) / prandtl(ig))**(1./2.)),f_ri_cd_min) |
---|
252 | fh(ig) = max((fm(ig) / prandtl(ig)), f_ri_cd_min) |
---|
253 | ELSE |
---|
254 | ! UNSTABLE BOUNDARY LAYER : |
---|
255 | sm = 2./rpi * (cinf - cn) * atan(-rib(ig)/ri0) + cn |
---|
256 | prandtl(ig) = -2./rpi * (pr_asym - pr_neut) * atan(rib(ig)/ri1) + pr_neut |
---|
257 | fm(ig) = MAX((sm**(3./2.) * sqrt(cepsilon) * (1 - rib(ig) / prandtl(ig))**(1./2.)),f_ri_cd_min) |
---|
258 | fh(ig) = MAX((fm(ig) / prandtl(ig)), f_ri_cd_min) |
---|
259 | ENDIF ! Rib < 0 |
---|
260 | ELSE |
---|
261 | ! Computation following parameterizaiton from from D.E. England et al. (95) |
---|
262 | IF (rib(ig) .gt. 0.) THEN |
---|
263 | ! STABLE BOUNDARY LAYER : |
---|
264 | prandtl(ig) = 1. |
---|
265 | IF(rib(ig) .lt. ric_colaitis) THEN |
---|
266 | ! Assuming alphah=1. and bh=bm for stable conditions : |
---|
267 | fm(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2 |
---|
268 | fh(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2 |
---|
269 | ELSE |
---|
270 | ! For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface |
---|
271 | fm(ig) = 1. |
---|
272 | fh(ig) = 1. |
---|
273 | ENDIF |
---|
274 | ELSE |
---|
275 | ! UNSTABLE BOUNDARY LAYER : |
---|
276 | fm(ig) = sqrt(1.-lambda*bm*rib(ig)) |
---|
277 | fh(ig) = (1./alphah)*((1.-lambda*bh*rib(ig))**0.5)*(1.-lambda*bm*rib(ig))**0.25 |
---|
278 | prandtl(ig) = alphah*((1.-lambda*bm*rib(ig))**0.25)/((1.-lambda*bh*rib(ig))**0.5) |
---|
279 | ENDIF ! Rib < 0 |
---|
280 | ENDIF ! atke |
---|
281 | ! Recompute the reynolds and z0t |
---|
282 | reynolds(ig) = karman*sqrt(fm(ig))*sqrt(zu2(ig))*pz0(ig)/(log(z1z0)*nu) |
---|
283 | pz0tcomp(ig) = pz0(ig)*exp(-karman*7.3*(reynolds(ig)**0.25)*(prandtl(ig)**0.5)+5*karman) |
---|
284 | residual = abs(pz0t-pz0tcomp(ig)) |
---|
285 | ite = ite+1 |
---|
286 | ENDDO ! of while |
---|
287 | ! Compute the coefficients Cdv; Cdh : neutral coefficient x stability functions |
---|
288 | pcdv(ig) = cdn(ig)*fm(ig) |
---|
289 | pcdh(ig) = chn(ig)*fh(ig) |
---|
290 | pz0t = 0. ! for next grid point |
---|
291 | ENDDO ! of ngrid |
---|
292 | |
---|
293 | !------------------------------------------------------------------------ |
---|
294 | !------------------------------------------------------------------------ |
---|
295 | ! PART II : USTAR/TSTAR/U_OUT/TETA_OUT COMPUTATIONS |
---|
296 | !------------------------------------------------------------------------ |
---|
297 | !------------------------------------------------------------------------ |
---|
298 | |
---|
299 | ! u* theta* computation |
---|
300 | DO n=1,n_out |
---|
301 | zout = z_out(n) |
---|
302 | DO ig=1,ngrid |
---|
303 | IF (rib(ig) .ge. ric_4interp) THEN !stable case, no turbulence |
---|
304 | ustar(ig) = 0. |
---|
305 | tstar(ig) = 0. |
---|
306 | ELSE |
---|
307 | ustar(ig) = sqrt(pcdv(ig))*sqrt(zu2(ig)) ! By definition, u* = sqrt(tau/U^2) = sqrt(Cdv) |
---|
308 | tstar(ig) = -pcdh(ig)*(pts(ig)-ph(ig,1))/sqrt(pcdv(ig)) ! theta* = (T1-Tsurf)*Cdh(ig)/u* |
---|
309 | ENDIF |
---|
310 | |
---|
311 | ! Interpolation: |
---|
312 | |
---|
313 | IF(zout .lt. pz0tcomp(ig)) THEN ! If z_out is in the surface layer , theta = tsurf; u = 0 |
---|
314 | u_out(ig,n) = 0. |
---|
315 | Teta_out(ig,n) = pts(ig) |
---|
316 | ELSEIF (zout .lt. pz0(ig)) THEN |
---|
317 | u_out(ig,n) = 0. |
---|
318 | ELSE |
---|
319 | IF (rib(ig) .ge. ric_4interp) THEN ! ustar=tstar=0 (and fm=fh=0) because no turbulence |
---|
320 | u_out(ig,n) = 0 |
---|
321 | Teta_out(ig,n) = pts(ig) |
---|
322 | ELSE |
---|
323 | ! We use the MO theory: du/dz = u*/(kappa z) *1 /sqrt(fm); dtheta/zd = theta*/(Pr kappa z) * fm/fh |
---|
324 | u_out(ig,n)= ustar(ig)*log(zout/pz0(ig))/(karman*sqrt(fm(ig))) |
---|
325 | Teta_out(ig,n)=pts(ig)+(tstar(ig)*sqrt(fm(ig))*log(zout/(pz0tcomp(ig)))/(karman*fh(ig))) |
---|
326 | ENDIF |
---|
327 | ENDIF |
---|
328 | ENDDO ! loop on n_out |
---|
329 | |
---|
330 | ! when using convective adjustment without thermals, a vertical potential temperature |
---|
331 | ! profile is assumed up to the thermal roughness length. Hence, theoretically, theta |
---|
332 | ! interpolated at any height in the surface layer is theta at the first level. |
---|
333 | |
---|
334 | IF ((.not.calltherm).and.(calladj)) THEN |
---|
335 | Teta_out(:,n)=ph(:,1) |
---|
336 | u_out(:,n)=(sqrt(cdn(:))*sqrt(zu2(ig))/karman)*log(zout/pz0(:)) |
---|
337 | ENDIF |
---|
338 | T_out(:,n) = Teta_out(:,n)*(exp((zout/zzlay(:,1))*(log(pplay(:,1)/ps))))**rcp ! not sure of what is done here: hydrostatic correction to account for the difference of pressure between zout and z1 ? |
---|
339 | ENDDO !of n=1,n_out |
---|
340 | |
---|
341 | |
---|
342 | !------------------------------------------------------------------------ |
---|
343 | !------------------------------------------------------------------------ |
---|
344 | ! PART III : VERTICAL_VELOCITY_VARIANCE/VERTICAL_TURBULENT_FLUX PROFILES |
---|
345 | !------------------------------------------------------------------------ |
---|
346 | !------------------------------------------------------------------------ |
---|
347 | |
---|
348 | ! We follow Spiga et. al 2010 (QJRMS) |
---|
349 | ! ------------ |
---|
350 | IF (calltherm) THEN |
---|
351 | DO ig=1, ngrid |
---|
352 | IF (zmax(ig) .gt. 0.) THEN |
---|
353 | x(ig) = zout/zmax(ig) |
---|
354 | ELSE |
---|
355 | x(ig) = 999. |
---|
356 | ENDIF |
---|
357 | ENDDO |
---|
358 | |
---|
359 | DO ig=1, ngrid |
---|
360 | ! dimensionless vertical heat flux |
---|
361 | IF (x(ig) .le. 0.3) THEN |
---|
362 | dvhf(ig) = ((-3.85/log(x(ig)))+0.07*log(x(ig))) *exp(-4.61*x(ig)) |
---|
363 | ELSEIF (x(ig) .le. 1.) THEN |
---|
364 | dvhf(ig) = -1.52*x(ig) + 1.24 |
---|
365 | ELSE |
---|
366 | dvhf(ig) = 0. |
---|
367 | ENDIF |
---|
368 | ! dimensionless vertical velocity variance |
---|
369 | IF (x(ig) .le. 1.) THEN |
---|
370 | dvvv(ig) = 2.05*(x(ig)**(2./3.))*(1.-0.64*x(ig))**2 |
---|
371 | ELSE |
---|
372 | dvvv(ig) = 0. |
---|
373 | ENDIF |
---|
374 | ENDDO |
---|
375 | |
---|
376 | vhf(:) = dvhf(:)*hfmax(:) |
---|
377 | vvv(:) = dvvv(:)*(wstar_in(:))**2 |
---|
378 | |
---|
379 | ENDIF ! of if calltherm |
---|
380 | |
---|
381 | |
---|
382 | !------------------------------------------------------------------------ |
---|
383 | !------------------------------------------------------------------------ |
---|
384 | ! PART IV : Outputs |
---|
385 | !------------------------------------------------------------------------ |
---|
386 | !------------------------------------------------------------------------ |
---|
387 | |
---|
388 | |
---|
389 | #ifndef MESOSCALE |
---|
390 | call write_output('tke_pbl','Tke in the pbl after physic','J/kg',tke(:,:)) |
---|
391 | call write_output('rib_pbl','Richardson in pbl parameter','m/s',rib(:)) |
---|
392 | call write_output('cdn_pbl','neutral momentum coef','m/s',cdn(:)) |
---|
393 | call write_output('fm_pbl','momentum stab function','m/s',fm(:)) |
---|
394 | call write_output('uv','wind norm first layer','m/s',sqrt(zu2(:))) |
---|
395 | call write_output('uvtrue','wind norm first layer','m/s',sqrt(pu(:,1)**2+pv(:,1)**2)) |
---|
396 | call write_output('chn_pbl','neutral momentum coef','m/s',chn(:)) |
---|
397 | call write_output('fh_pbl','momentum stab function','m/s',fh(:)) |
---|
398 | call write_output('B_pbl','buoyancy','m/',(ph(:,1)-pts(:))/pts(:)) |
---|
399 | call write_output('zot_pbl','buoyancy','ms',pz0tcomp(:)) |
---|
400 | call write_output('zz1','1st layer altitude','m',zzlay(:,1)) |
---|
401 | #endif |
---|
402 | |
---|
403 | RETURN |
---|
404 | END SUBROUTINE pbl_parameters |
---|
405 | |
---|
406 | END MODULE pbl_parameters_mod |
---|
407 | |
---|