Line data Source code
1 : !-------------------------------------------------------------------------
2 : ! $Id$
3 : !===============================================================================
4 : module sfc_varnce_module
5 :
6 : implicit none
7 :
8 : private ! Default to private
9 :
10 : public :: calc_sfc_varnce
11 :
12 : contains
13 :
14 : !=============================================================================
15 352944 : subroutine calc_sfc_varnce( nz, ngrdcol, gr, dt, sfc_elevation, &
16 352944 : upwp_sfc, vpwp_sfc, wpthlp, wprtp_sfc, &
17 352944 : um, vm, Lscale_up, wpsclrp_sfc, &
18 352944 : lhs_splat_wp2, tau_zm, &
19 : !wp2_splat_sfc, tau_zm_sfc, &
20 : l_vary_convect_depth, &
21 : clubb_params, &
22 : stats_metadata, &
23 352944 : stats_zm, &
24 352944 : wp2, up2, vp2, &
25 352944 : thlp2, rtp2, rtpthlp, &
26 352944 : sclrp2, sclrprtp, sclrpthlp )
27 :
28 : ! Description:
29 : ! This subroutine computes estimate of the surface thermodynamic and wind
30 : ! component second order moments.
31 :
32 : ! References:
33 : ! Andre, J. C., G. De Moor, P. Lacarrere, G. Therry, and R. Du Vachat, 1978.
34 : ! Modeling the 24-Hour Evolution of the Mean and Turbulent Structures of
35 : ! the Planetary Boundary Layer. J. Atmos. Sci., 35, 1861 -- 1883.
36 : !-----------------------------------------------------------------------
37 :
38 : use grid_class, only: &
39 : grid
40 :
41 : use parameters_model, only: &
42 : T0 ! Variable(s)
43 :
44 : use constants_clubb, only: &
45 : four, & ! Variable(s)
46 : two, &
47 : one, &
48 : two_thirds, &
49 : one_third, &
50 : one_fourth, &
51 : zero, &
52 : grav, &
53 : eps, &
54 : w_tol_sqd, &
55 : thl_tol, &
56 : rt_tol, &
57 : max_mag_correlation_flux, &
58 : fstderr, &
59 : wp2_max
60 :
61 : use parameters_model, only: &
62 : sclr_dim ! Variable(s)
63 :
64 : use parameter_indices, only: &
65 : nparams, &
66 : ia_const, &
67 : iup2_sfc_coef
68 :
69 : use numerical_check, only: &
70 : sfc_varnce_check ! Procedure
71 :
72 : use error_code, only: &
73 : clubb_at_least_debug_level, & ! Procedure
74 : err_code, & ! Error Indicator
75 : clubb_fatal_error ! Constant
76 :
77 : use array_index, only: &
78 : iisclr_rt, & ! Index for a scalar emulating rt
79 : iisclr_thl ! Index for a scalar emulating thetal
80 :
81 : use stats_type, only: &
82 : stats ! Type
83 :
84 : use stats_type_utilities, only: &
85 : stat_end_update_pt, & ! Procedure(s)
86 : stat_begin_update_pt, &
87 : stat_update_var_pt
88 :
89 : use stats_variables, only: &
90 : stats_metadata_type
91 :
92 : use clubb_precision, only: &
93 : core_rknd ! Variable(s)
94 :
95 : implicit none
96 :
97 : ! Constant Parameters
98 :
99 : ! Logical for Andre et al., 1978 parameterization.
100 : logical, parameter :: l_andre_1978 = .false.
101 :
102 : real( kind = core_rknd ), parameter :: &
103 : z_const = one, & ! Defined height of 1 meter [m]
104 : ! Vince Larson increased ufmin to stabilize arm_97. 24 Jul 2007
105 : ! ufmin = 0.0001_core_rknd, &
106 : ufmin = 0.01_core_rknd, & ! Minimum allowable value of u* [m/s]
107 : ! End Vince Larson's change.
108 : ! Vince Larson changed in order to make correlations between [-1,1]. 31 Jan 2008.
109 : ! sclr_var_coef = 0.25_core_rknd, & ! This value is made up! - Vince Larson 12 Jul 2005
110 : sclr_var_coef = 0.4_core_rknd, & ! This value is made up! - Vince Larson 12 Jul 2005
111 : ! End Vince Larson's change
112 : ! Vince Larson reduced surface spike in scalar variances associated
113 : ! w/ Andre et al. 1978 scheme
114 : reduce_coef = 0.2_core_rknd
115 :
116 : !-------------------------- Input Variables --------------------------
117 : integer, intent(in) :: &
118 : nz, &
119 : ngrdcol
120 :
121 : type(grid), intent(in) :: &
122 : gr
123 :
124 : real( kind = core_rknd ) :: &
125 : dt
126 :
127 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
128 : sfc_elevation, &
129 : upwp_sfc, & ! Surface u momentum flux, <u'w'>|_sfc [m^2/s^2]
130 : vpwp_sfc, & ! Surface v momentum flux, <v'w'>|_sfc [m^2/s^2]
131 : wprtp_sfc ! Surface moisture flux, <w'rt'>|_sfc [kg/kg m/s]
132 :
133 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
134 : wpthlp, & ! Surface thetal flux, <w'thl'>|_sfc [K m/s]
135 : um, & ! Surface u wind component, <u> [m/s]
136 : vm, & ! Surface v wind component, <v> [m/s]
137 : Lscale_up, & ! Upward component of Lscale at surface [m]
138 : lhs_splat_wp2, &
139 : !wp2_splat, & ! Tendency of <w'^2> due to splatting of eddies at zm(1) [m^2/s^3]
140 : tau_zm ! Turbulent dissipation time at level zm(1) [s]
141 :
142 :
143 : real( kind = core_rknd ), dimension(ngrdcol,sclr_dim), intent(in) :: &
144 : wpsclrp_sfc ! Passive scalar flux, <w'sclr'>|_sfc [units m/s]
145 :
146 : logical, intent(in) :: &
147 : l_vary_convect_depth
148 :
149 : real( kind = core_rknd ), dimension(nparams), intent(in) :: &
150 : clubb_params ! Array of CLUBB's tunable parameters [units vary]
151 :
152 : type (stats_metadata_type), intent(in) :: &
153 : stats_metadata
154 :
155 : !-------------------------- InOut Variables --------------------------
156 : type (stats), target, intent(inout), dimension(ngrdcol) :: &
157 : stats_zm
158 :
159 : !-------------------------- Output Variables --------------------------
160 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(inout) :: &
161 : wp2, & ! Surface variance of w, <w'^2>|_sfc [m^2/s^2]
162 : up2, & ! Surface variance of u, <u'^2>|_sfc [m^2/s^2]
163 : vp2, & ! Surface variance of v, <v'^2>|_sfc [m^2/s^2]
164 : thlp2, & ! Surface variance of theta-l, <thl'^2>|_sfc [K^2]
165 : rtp2, & ! Surface variance of rt, <rt'^2>|_sfc [(kg/kg)^2]
166 : rtpthlp ! Surface covariance of rt and theta-l [kg K/kg]
167 :
168 : real( kind = core_rknd ), dimension(ngrdcol,nz,sclr_dim), intent(inout) :: &
169 : sclrp2, & ! Surface variance of passive scalar [units^2]
170 : sclrprtp, & ! Surface covariance of pssv scalar and rt [units kg/kg]
171 : sclrpthlp ! Surface covariance of pssv scalar and theta-l [units K]
172 :
173 : !-------------------------- Local Variables --------------------------
174 : real( kind = core_rknd ), dimension(ngrdcol) :: &
175 705888 : uf, & ! Surface friction vel., u*, in older version [m/s]
176 705888 : depth_pos_wpthlp, & ! Thickness of the layer near the surface with wpthlp > 0 [m]
177 705888 : min_wp2_sfc_val ! Minimum value of wp2_sfc that guarantees [m^2/s^2]
178 : ! correlation of (w,rt) and (w,thl) is within (-1,1)
179 :
180 : ! Variables for Andre et al., 1978 parameterization.
181 : real( kind = core_rknd ), dimension(ngrdcol) :: &
182 352944 : um_sfc_sqd, & ! Surface value of <u>^2 [m^2/s^2]
183 352944 : vm_sfc_sqd, & ! Surface value of <v>^2 [m^2/s^2]
184 352944 : usp2_sfc, & ! u_s (vector oriented w/ mean sfc. wind) variance [m^2/s^2]
185 352944 : vsp2_sfc, & ! v_s (vector perpen. to mean sfc. wind) variance [m^2/s^2]
186 352944 : ustar, & ! Surface friction velocity, u* [m/s]
187 352944 : zeta, & ! Dimensionless height z_const/Lngth, where z_const = 1 m. [-]
188 705888 : wp2_splat_sfc_correction ! Reduction in wp2_sfc due to splatting [m^2/s^2]
189 :
190 : real( kind = core_rknd ) :: &
191 : ustar2, & ! Square of surface friction velocity, u* [m^2/s^2]
192 : wstar ! Convective velocity, w* [m/s]
193 :
194 : real( kind = core_rknd ) :: &
195 : Lngth ! Monin-Obukhov length [m]
196 :
197 : real( kind = core_rknd ) :: &
198 : up2_sfc_coef, & ! CLUBB tunable parameter up2_sfc_coef [-]
199 : a_const ! Coefficient in front of wp2, up2, and vp2
200 :
201 : integer :: i, k, sclr ! Loop index
202 :
203 : !-------------------------- Begin Code --------------------------
204 :
205 : !$acc enter data create( uf, depth_pos_wpthlp, min_wp2_sfc_val, &
206 : !$acc um_sfc_sqd, vm_sfc_sqd, usp2_sfc, vsp2_sfc, &
207 : !$acc ustar, zeta, wp2_splat_sfc_correction )
208 :
209 352944 : up2_sfc_coef = clubb_params(iup2_sfc_coef)
210 :
211 : ! Reflect surface varnce changes in budget
212 352944 : if ( stats_metadata%l_stats_samp ) then
213 :
214 : !$acc update host( wp2, up2, vp2, thlp2, rtp2, rtpthlp )
215 :
216 0 : do i = 1, ngrdcol
217 : call stat_begin_update_pt( stats_metadata%ithlp2_sf, 1, & ! intent(in)
218 0 : thlp2(i,1) / dt, & ! intent(in)
219 0 : stats_zm(i) ) ! intent(inout)
220 :
221 : call stat_begin_update_pt( stats_metadata%irtp2_sf, 1, & ! intent(in)
222 : rtp2(i,1) / dt, & ! intent(in)
223 0 : stats_zm(i) ) ! intent(inout)
224 :
225 : call stat_begin_update_pt( stats_metadata%irtpthlp_sf, 1, & ! intent(in)
226 : rtpthlp(i,1) / dt, & ! intent(in)
227 0 : stats_zm(i) ) ! intent(inout)
228 :
229 : call stat_begin_update_pt( stats_metadata%iup2_sf, 1, & ! intent(in)
230 : up2(i,1) / dt, & ! intent(in)
231 0 : stats_zm(i) ) ! intent(inout)
232 :
233 : call stat_begin_update_pt( stats_metadata%ivp2_sf, 1, & ! intent(in)
234 : vp2(i,1) / dt, & ! intent(in)
235 0 : stats_zm(i) ) ! intent(inout)
236 :
237 : call stat_begin_update_pt( stats_metadata%iwp2_sf, 1, & ! intent(in)
238 : wp2(i,1) / dt, & ! intent(in)
239 0 : stats_zm(i) ) ! intent(inout)
240 : end do
241 : end if
242 :
243 : !$acc parallel loop gang vector default(present)
244 5893344 : do i = 1, ngrdcol
245 : ! Find thickness of layer near surface with positive heat flux.
246 : ! This is used when l_vary_convect_depth=.true. in order to determine wp2.
247 5893344 : if ( wpthlp(i,1) <= zero ) then
248 1627158 : depth_pos_wpthlp(i) = one ! When sfc heat flux is negative, set depth to 1 m.
249 : else ! When sfc heat flux is positive, march up sounding until wpthlp 1st becomes negative.
250 : k = 1
251 23799670 : do while ( wpthlp(i,k) > zero .and. (gr%zm(i,k)-sfc_elevation(i)) < 1000._core_rknd )
252 19886428 : k = k + 1
253 : end do
254 3913242 : depth_pos_wpthlp(i) = max( one, gr%zm(i,k)-sfc_elevation(i) )
255 : end if
256 : end do
257 : !$acc end parallel loop
258 :
259 : ! a_const used to be set here; now it is a tunable parameter
260 : !if ( .not. l_vary_convect_depth ) then
261 : ! a_const = 1.8_core_rknd
262 : !else
263 : ! a_const = 0.6_core_rknd
264 : !end if
265 :
266 352944 : a_const = clubb_params(ia_const)
267 :
268 : if ( l_andre_1978 ) then
269 :
270 : ! Calculate <u>^2 and <v>^2.
271 : !$acc parallel loop gang vector default(present)
272 : do i = 1, ngrdcol
273 : um_sfc_sqd(i) = um(i,2)**2
274 : vm_sfc_sqd(i) = vm(i,2)**2
275 : end do
276 : !$acc end parallel loop
277 :
278 : ! Calculate surface friction velocity, u*.
279 : !$acc parallel loop gang vector default(present)
280 : do i = 1, ngrdcol
281 : ustar(i) = max( ( upwp_sfc(i)**2 + vpwp_sfc(i)**2 )**(one_fourth), ufmin )
282 : end do
283 : !$acc end parallel loop
284 :
285 : !$acc parallel loop gang vector default(present)
286 : do i = 1, ngrdcol
287 : if ( abs(wpthlp(i,1)) > eps) then
288 :
289 : ! Find Monin-Obukhov Length (Andre et al., 1978, p. 1866).
290 : Lngth = - ( ustar(i)**3 ) &
291 : / ( 0.35_core_rknd * (one/T0) * grav * wpthlp(i,1) )
292 :
293 : ! Find the value of dimensionless height zeta
294 : ! (Andre et al., 1978, p. 1866).
295 : zeta(i) = z_const / Lngth
296 :
297 : else ! wpthlp = 0
298 :
299 : ! The value of Monin-Obukhov length is +inf when ustar < 0 and -inf
300 : ! when ustar > 0. Either way, zeta = 0.
301 : zeta(i) = zero
302 :
303 : endif ! wpthlp /= 0
304 : end do
305 : !$acc end parallel loop
306 :
307 : ! Andre et al, 1978, Eq. 29.
308 : ! Notes: 1) "reduce_coef" is a reduction coefficient intended to make
309 : ! the values of rtp2, thlp2, and rtpthlp smaller at the
310 : ! surface.
311 : ! 2) With the reduction coefficient having a value of 0.2, the
312 : ! surface correlations of both w & rt and w & thl have a value
313 : ! of about 0.845. These correlations are greater if zeta < 0.
314 : ! The correlations have a value greater than 1 if
315 : ! zeta <= -0.212.
316 : ! 3) The surface correlation of rt & thl is 1.
317 : ! Brian Griffin; February 2, 2008.
318 : !$acc parallel loop gang vector default(present)
319 : do i = 1, ngrdcol
320 : if ( zeta(i) < zero ) then
321 :
322 : thlp2(i,1) = reduce_coef &
323 : * ( wpthlp(i,1)**2 / ustar(i)**2 ) &
324 : * four * ( one - 8.3_core_rknd * zeta(i) )**(-two_thirds)
325 :
326 : rtp2(i,1) = reduce_coef &
327 : * ( wprtp_sfc(i)**2 / ustar(i)**2 ) &
328 : * four * ( one - 8.3_core_rknd * zeta(i) )**(-two_thirds)
329 :
330 : rtpthlp(i,1) = reduce_coef &
331 : * ( wprtp_sfc(i) * wpthlp(i,1) / ustar(i)**2 ) &
332 : * four * ( one - 8.3_core_rknd * zeta(i) )**(-two_thirds)
333 :
334 : wp2(i,1) = ( ustar(i)**2 ) &
335 : * ( 1.75_core_rknd + two * (-zeta(i))**(two_thirds) )
336 :
337 : else
338 :
339 : thlp2(i,1) = reduce_coef &
340 : * four * ( wpthlp(i,1)**2 / ustar(i)**2 )
341 :
342 : rtp2(i,1) = reduce_coef &
343 : * four * ( wprtp_sfc(i)**2 / ustar(i)**2 )
344 :
345 : rtpthlp(i,1) = reduce_coef &
346 : * four * ( wprtp_sfc(i) * wpthlp(i,1) / ustar(i)**2 )
347 :
348 : wp2(i,1) = 1.75_core_rknd * ustar(i)**2
349 :
350 : endif
351 : end do
352 : !$acc end parallel loop
353 :
354 : !$acc parallel loop gang vector default(present)
355 : do i = 1, ngrdcol
356 : thlp2(i,1) = max( thl_tol**2, thlp2(i,1) )
357 : rtp2(i,1) = max( rt_tol**2, rtp2(i,1) )
358 : end do
359 : !$acc end parallel loop
360 :
361 : ! Calculate wstar following Andre et al., 1978, p. 1866.
362 : ! w* = ( ( 1 / T0 ) * g * <w'thl'>|_sfc * z_i )^(1/3);
363 : ! where z_i is the height of the mixed layer. The value of CLUBB's
364 : ! upward component of mixing length, Lscale_up, at the surface will be
365 : ! used as z_i.
366 : !$acc parallel loop gang vector default(present)
367 : do i = 1, ngrdcol
368 : wstar = ( (one/T0) * grav * wpthlp(i,1) * Lscale_up(i,2) )**(one_third)
369 :
370 : ! Andre et al., 1978, Eq. 29.
371 : ! Andre et al. (1978) defines horizontal wind surface variances in terms
372 : ! of orientation with the mean surface wind. The vector u_s is the wind
373 : ! vector oriented with the mean surface wind. The vector v_s is the wind
374 : ! vector oriented perpendicular to the mean surface wind. Thus, <u_s> is
375 : ! equal to the mean surface wind (both in speed and direction), and <v_s>
376 : ! is 0. Equation 29 gives the formula for the variance of u_s, which is
377 : ! <u_s'^2> (usp2_sfc in the code), and the formula for the variance of
378 : ! v_s, which is <v_s'^2> (vsp2_sfc in the code).
379 : if ( wpthlp(i,1) > zero ) then
380 :
381 : usp2_sfc(i) = four * ustar(i)**2 + 0.3_core_rknd * wstar**2
382 :
383 : vsp2_sfc(i) = 1.75_core_rknd * ustar(i)**2 + 0.3_core_rknd * wstar**2
384 :
385 : else
386 :
387 : usp2_sfc(i) = four * ustar(i)**2
388 :
389 : vsp2_sfc(i) = 1.75_core_rknd * ustar(i)**2
390 :
391 : endif
392 : end do
393 : !$acc end parallel loop
394 :
395 : ! Add effect of vertical compression of eddies on horizontal gustiness.
396 : ! First, ensure that wp2 does not make the correlation
397 : ! of (w,thl) or (w,rt) outside (-1,1).
398 : ! Perhaps in the future we should also ensure that the correlations
399 : ! of (w,u) and (w,v) are not outside (-1,1).
400 : !$acc parallel loop gang vector default(present)
401 : do i = 1, ngrdcol
402 : min_wp2_sfc_val(i) &
403 : = max( w_tol_sqd, &
404 : wprtp_sfc(i)**2 / ( rtp2(i,1) * max_mag_correlation_flux**2 ), &
405 : wpthlp(i,1)**2 / ( thlp2(i,1) * max_mag_correlation_flux**2 ) )
406 : end do
407 : !$acc end parallel loop
408 :
409 : !$acc parallel loop gang vector default(present)
410 : do i = 1, ngrdcol
411 : if ( wp2(i,1) - tau_zm(i,1) * lhs_splat_wp2(i,1) * wp2(i,1) &
412 : < min_wp2_sfc_val(i) ) then
413 : !if ( wp2(i,1) + tau_zm(i,1) * wp2_splat(i,1) < min_wp2_sfc_val(i) ) then
414 : ! splatting correction drives wp2 to overly small value
415 : wp2_splat_sfc_correction(i) = -wp2(i,1) + min_wp2_sfc_val(i)
416 : wp2(i,1) = min_wp2_sfc_val(i)
417 : else
418 : wp2_splat_sfc_correction(i) = tau_zm(i,1) * lhs_splat_wp2(i,1) * wp2(i,1)
419 : !wp2_splat_sfc_correction(i) = tau_zm(i,1) * wp2_splat(i,1)
420 : wp2(i,1) = wp2(i,1) + wp2_splat_sfc_correction(i)
421 : end if
422 : end do
423 : !$acc end parallel loop
424 :
425 : !$acc parallel loop gang vector default(present)
426 : do i = 1, ngrdcol
427 : usp2_sfc(i) = usp2_sfc(i) - 0.5_core_rknd * wp2_splat_sfc_correction(i)
428 : vsp2_sfc(i) = vsp2_sfc(i) - 0.5_core_rknd * wp2_splat_sfc_correction(i)
429 : end do
430 : !$acc end parallel loop
431 :
432 : ! Variance of u, <u'^2>, at the surface can be found from <u_s'^2>,
433 : ! <v_s'^2>, and mean winds (at the surface) <u> and <v>, such that:
434 : ! <u'^2>|_sfc = <u_s'^2> * [ <u>^2 / ( <u>^2 + <v>^2 ) ]
435 : ! + <v_s'^2> * [ <v>^2 / ( <u>^2 + <v>^2 ) ];
436 : ! where <u>^2 + <v>^2 /= 0.
437 : !$acc parallel loop gang vector default(present)
438 : do i = 1, ngrdcol
439 : up2(i,1) = usp2_sfc(i) * ( um_sfc_sqd(i) / max( um_sfc_sqd(i) + vm_sfc_sqd(i) , eps ) ) &
440 : + vsp2_sfc(i) * ( vm_sfc_sqd(i) / max( um_sfc_sqd(i) + vm_sfc_sqd(i) , eps ) )
441 : end do
442 : !$acc end parallel loop
443 :
444 : ! Variance of v, <v'^2>, at the surface can be found from <u_s'^2>,
445 : ! <v_s'^2>, and mean winds (at the surface) <u> and <v>, such that:
446 : ! <v'^2>|_sfc = <v_s'^2> * [ <u>^2 / ( <u>^2 + <v>^2 ) ]
447 : ! + <u_s'^2> * [ <v>^2 / ( <u>^2 + <v>^2 ) ];
448 : ! where <u>^2 + <v>^2 /= 0.
449 : !$acc parallel loop gang vector default(present)
450 : do i = 1, ngrdcol
451 : vp2(i,1) = vsp2_sfc(i) * ( um_sfc_sqd(i) / max( um_sfc_sqd(i) + vm_sfc_sqd(i) , eps ) ) &
452 : + usp2_sfc(i) * ( vm_sfc_sqd(i) / max( um_sfc_sqd(i) + vm_sfc_sqd(i) , eps ) )
453 : end do
454 : !$acc end parallel loop
455 :
456 : ! Passive scalars
457 : if ( sclr_dim > 0 ) then
458 :
459 : !$acc parallel loop gang vector collapse(2) default(present)
460 : do sclr = 1, sclr_dim
461 : do i = 1, ngrdcol
462 :
463 : ! Notes: 1) "reduce_coef" is a reduction coefficient intended to
464 : ! make the values of sclrprtp, sclrpthlp, and sclrp2
465 : ! smaller at the surface.
466 : ! 2) With the reduction coefficient having a value of 0.2,
467 : ! the surface correlation of w & sclr has a value of
468 : ! about 0.845. The correlation is greater if zeta < 0.
469 : ! The correlation has a value greater than 1 if
470 : ! zeta <= -0.212.
471 : ! 3) The surface correlations of both rt & sclr and
472 : ! thl & sclr are 1.
473 : ! Brian Griffin; February 2, 2008.
474 : if ( zeta(i) < zero ) then
475 :
476 : sclrprtp(i,1,sclr) &
477 : = reduce_coef &
478 : * ( wpsclrp_sfc(i,sclr) * wprtp_sfc(i) / ustar(i)**2 ) &
479 : * four * ( one - 8.3_core_rknd * zeta(i) )**(-two_thirds)
480 :
481 : sclrpthlp(i,1,sclr) &
482 : = reduce_coef &
483 : * ( wpsclrp_sfc(i,sclr) * wpthlp(i,1) / ustar(i)**2 ) &
484 : * four * ( one - 8.3_core_rknd * zeta(i) )**(-two_thirds)
485 :
486 : sclrp2(i,1,sclr) &
487 : = reduce_coef &
488 : * ( wpsclrp_sfc(i,sclr)**2 / ustar(i)**2 ) &
489 : * four * ( one - 8.3_core_rknd * zeta(i) )**(-two_thirds)
490 :
491 : else
492 :
493 : sclrprtp(i,1,sclr) &
494 : = reduce_coef &
495 : * four * ( wpsclrp_sfc(i,sclr) * wprtp_sfc(i) / ustar(i)**2 )
496 :
497 : sclrpthlp(i,1,sclr) &
498 : = reduce_coef &
499 : * four * ( wpsclrp_sfc(i,sclr) * wpthlp(i,1) / ustar(i)**2 )
500 :
501 : sclrp2(i,1,sclr) &
502 : = reduce_coef &
503 : * four * ( wpsclrp_sfc(i,sclr)**2 / ustar(i)**2 )
504 :
505 : endif
506 :
507 : end do
508 : enddo ! i = 1, sclr_dim
509 : !$acc end parallel loop
510 :
511 : endif
512 :
513 :
514 : else ! Previous code.
515 :
516 : ! Compute ustar^2
517 : !$acc parallel loop gang vector default(present)
518 5893344 : do i = 1, ngrdcol
519 5540400 : ustar2 = sqrt( upwp_sfc(i) * upwp_sfc(i) + vpwp_sfc(i) * vpwp_sfc(i) )
520 :
521 : ! Compute wstar following Andre et al., 1976
522 5540400 : if ( wpthlp(i,1) > zero ) then
523 3913242 : if ( .not. l_vary_convect_depth ) then
524 3913242 : wstar = ( one/T0 * grav * wpthlp(i,1) * z_const )**(one_third)
525 : else
526 0 : wstar = ( one/T0 * grav * wpthlp(i,1) * 0.2_core_rknd * depth_pos_wpthlp(i) )**(one_third)
527 : end if
528 : else
529 : wstar = zero
530 : endif
531 :
532 : ! Surface friction velocity following Andre et al. 1978
533 5540400 : if ( .not. l_vary_convect_depth ) then
534 5540400 : uf(i) = sqrt( ustar2 + 0.3_core_rknd * wstar * wstar )
535 : else
536 0 : uf(i) = sqrt( ustar2 + wstar * wstar )
537 : end if
538 :
539 5893344 : uf(i) = max( ufmin, uf(i) )
540 : end do
541 : !$acc end parallel loop
542 :
543 : ! Compute estimate for surface second order moments
544 : !$acc parallel loop gang vector default(present)
545 5893344 : do i = 1, ngrdcol
546 5540400 : wp2(i,1) = a_const * uf(i)**2
547 5540400 : up2(i,1) = up2_sfc_coef * a_const * uf(i)**2 ! From Andre, et al. 1978
548 5893344 : vp2(i,1) = up2_sfc_coef * a_const * uf(i)**2 ! " "
549 : end do
550 : !$acc end parallel loop
551 :
552 : ! Notes: 1) With "a" having a value of 1.8, the surface correlations of
553 : ! both w & rt and w & thl have a value of about 0.878.
554 : ! 2) The surface correlation of rt & thl is 0.5.
555 : ! Brian Griffin; February 2, 2008.
556 :
557 352944 : if ( .not. l_vary_convect_depth ) then
558 : !$acc parallel loop gang vector default(present)
559 5893344 : do i = 1, ngrdcol
560 5540400 : thlp2(i,1) = 0.4_core_rknd * a_const * ( wpthlp(i,1) / uf(i) )**2
561 5540400 : rtp2(i,1) = 0.4_core_rknd * a_const * ( wprtp_sfc(i) / uf(i) )**2
562 : rtpthlp(i,1) = 0.2_core_rknd * a_const * ( wpthlp(i,1) / uf(i) ) &
563 5893344 : * ( wprtp_sfc(i) / uf(i) )
564 : end do
565 : !$acc end parallel loop
566 : else
567 : !$acc parallel loop gang vector default(present)
568 0 : do i = 1, ngrdcol
569 0 : thlp2(i,1) = ( wpthlp(i,1) / uf(i) )**2 / ( max_mag_correlation_flux**2 * a_const )
570 0 : rtp2(i,1) = ( wprtp_sfc(i) / uf(i) )**2 / ( max_mag_correlation_flux**2 * a_const )
571 0 : rtpthlp(i,1) = max_mag_correlation_flux * sqrt( thlp2(i,1) * rtp2(i,1) )
572 : end do
573 : !$acc end parallel loop
574 : end if
575 :
576 : !$acc parallel loop gang vector default(present)
577 5893344 : do i = 1, ngrdcol
578 5540400 : thlp2(i,1) = max( thl_tol**2, thlp2(i,1) )
579 5893344 : rtp2(i,1) = max( rt_tol**2, rtp2(i,1) )
580 : end do
581 : !$acc end parallel loop
582 :
583 : ! Add effect of vertical compression of eddies on horizontal gustiness.
584 : ! First, ensure that wp2 does not make the correlation
585 : ! of (w,thl) or (w,rt) outside (-1,1).
586 : ! Perhaps in the future we should also ensure that the correlations
587 : ! of (w,u) and (w,v) are not outside (-1,1).
588 : !$acc parallel loop gang vector default(present)
589 5893344 : do i = 1, ngrdcol
590 5540400 : min_wp2_sfc_val(i) &
591 : = max( w_tol_sqd, &
592 0 : wprtp_sfc(i)**2 / ( rtp2(i,1) * max_mag_correlation_flux**2 ), &
593 11433744 : wpthlp(i,1)**2 / ( thlp2(i,1) * max_mag_correlation_flux**2 ) )
594 : end do
595 : !$acc end parallel loop
596 :
597 : !$acc parallel loop gang vector default(present)
598 5893344 : do i = 1, ngrdcol
599 5540400 : if ( wp2(i,1) - tau_zm(i,1) * lhs_splat_wp2(i,1) * wp2(i,1) &
600 : < min_wp2_sfc_val(i) ) then
601 : !if ( wp2(i,1) + tau_zm(i,1) * wp2_splat(i,1) < min_wp2_sfc_val(i) ) then
602 : ! splatting correction drives wp2 to overly small values
603 4132 : wp2_splat_sfc_correction(i) = -wp2(i,1) + min_wp2_sfc_val(i)
604 4132 : wp2(i,1) = min_wp2_sfc_val(i)
605 : else
606 5536268 : wp2_splat_sfc_correction(i) = tau_zm(i,1) * lhs_splat_wp2(i,1) * wp2(i,1)
607 : !wp2_splat_sfc_correction(i) = tau_zm(i,1) * wp2_splat(i,1)
608 5536268 : wp2(i,1) = wp2(i,1) + wp2_splat_sfc_correction(i)
609 : end if
610 :
611 5540400 : up2(i,1) = up2(i,1) - 0.5_core_rknd * wp2_splat_sfc_correction(i)
612 5893344 : vp2(i,1) = vp2(i,1) - 0.5_core_rknd * wp2_splat_sfc_correction(i)
613 : end do
614 : !$acc end parallel loop
615 :
616 : ! Passive scalars
617 352944 : if ( sclr_dim > 0 ) then
618 : !$acc parallel loop gang vector collapse(2) default(present)
619 0 : do sclr = 1, sclr_dim
620 0 : do i = 1, ngrdcol
621 : ! Vince Larson changed coeffs to make correlations between [-1,1].
622 : ! 31 Jan 2008
623 : ! sclrprtp(i) &
624 : ! = a * (wprtp_sfc / uf) * (wpsclrp(i) / uf)
625 : ! sclrpthlp(i) &
626 : ! = a * (wpthlp / uf) * (wpsclrp(i) / uf)
627 : ! sclrp2(i) &
628 : ! = sclr_var_coef * a * ( wpsclrp(i) / uf )**2
629 : ! Notes: 1) With "a" having a value of 1.8 and "sclr_var_coef"
630 : ! having a value of 0.4, the surface correlation of
631 : ! w & sclr has a value of about 0.878.
632 : ! 2) With "sclr_var_coef" having a value of 0.4, the
633 : ! surface correlations of both rt & sclr and
634 : ! thl & sclr are 0.5.
635 : ! Brian Griffin; February 2, 2008.
636 :
637 : ! We use the following if..then's to make sclr_rt and sclr_thl
638 : ! close to the actual thlp2/rtp2 at the surface.
639 : ! -dschanen 25 Sep 08
640 0 : if ( sclr == iisclr_rt ) then
641 : ! If we are trying to emulate rt with the scalar, then we
642 : ! use the variance coefficient from above
643 0 : sclrprtp(i,1,sclr) = 0.4_core_rknd * a_const &
644 0 : * ( wprtp_sfc(i) / uf(i) ) * ( wpsclrp_sfc(i,sclr) / uf(i) )
645 : else
646 0 : sclrprtp(i,1,sclr) = 0.2_core_rknd * a_const &
647 0 : * ( wprtp_sfc(i) / uf(i) ) * ( wpsclrp_sfc(i,sclr) / uf(i) )
648 : endif
649 :
650 0 : if ( sclr == iisclr_thl ) then
651 : ! As above, but for thetal
652 0 : sclrpthlp(i,1,sclr) = 0.4_core_rknd * a_const &
653 : * ( wpthlp(i,1) / uf(i) ) &
654 0 : * ( wpsclrp_sfc(i,sclr) / uf(i) )
655 : else
656 0 : sclrpthlp(i,1,sclr) = 0.2_core_rknd * a_const &
657 : * ( wpthlp(i,1) / uf(i) ) &
658 0 : * ( wpsclrp_sfc(i,sclr) / uf(i) )
659 : endif
660 :
661 0 : sclrp2(i,1,sclr) = sclr_var_coef * a_const &
662 0 : * ( wpsclrp_sfc(i,sclr) / uf(i) )**2
663 :
664 : ! End Vince Larson's change.
665 :
666 : end do
667 : enddo ! 1,...sclr_dim
668 : !$acc end parallel loop
669 :
670 : endif ! sclr_dim > 0
671 :
672 : endif ! l_andre_1978
673 :
674 : ! Clip wp2 at wp2_max, same as in advance_wp2_wp3
675 : !$acc parallel loop gang vector default(present)
676 5893344 : do i = 1, ngrdcol
677 5893344 : wp2(i,1) = min( wp2(i,1), wp2_max )
678 : end do
679 : !$acc end parallel loop
680 :
681 : !$acc parallel loop gang vector default(present)
682 5893344 : do i = 1, ngrdcol
683 5893344 : if ( abs(gr%zm(i,1)-sfc_elevation(i)) > abs(gr%zm(i,1)+sfc_elevation(i))*eps/2 ) then
684 :
685 : ! Variances for cases where the lowest level is not at the surface.
686 : ! Eliminate surface effects on lowest level variances.
687 0 : wp2(i,1) = w_tol_sqd
688 0 : up2(i,1) = w_tol_sqd
689 0 : vp2(i,1) = w_tol_sqd
690 0 : thlp2(i,1) = thl_tol**2
691 0 : rtp2(i,1) = rt_tol**2
692 0 : rtpthlp(i,1) = 0.0_core_rknd
693 :
694 : end if ! gr%zm(1,1) == sfc_elevation
695 : end do
696 : !$acc end parallel loop
697 :
698 352944 : if ( sclr_dim > 0 ) then
699 :
700 : !$acc parallel loop gang vector collapse(2) default(present)
701 0 : do sclr = 1, sclr_dim
702 0 : do i = 1, ngrdcol
703 0 : if ( abs(gr%zm(i,1)-sfc_elevation(i)) > abs(gr%zm(i,1)+sfc_elevation(i))*eps/2 ) then
704 :
705 : ! Variances for cases where the lowest level is not at the surface.
706 : ! Eliminate surface effects on lowest level variances.
707 0 : sclrp2(i,1,sclr) = 0.0_core_rknd
708 0 : sclrprtp(i,1,sclr) = 0.0_core_rknd
709 0 : sclrpthlp(i,1,sclr) = 0.0_core_rknd
710 : end if
711 : end do
712 : end do
713 : end if
714 :
715 352944 : if ( stats_metadata%l_stats_samp ) then
716 :
717 : !$acc update host( wp2, up2, vp2, thlp2, rtp2, rtpthlp )
718 :
719 0 : do i = 1, ngrdcol
720 : call stat_end_update_pt( stats_metadata%ithlp2_sf, 1, & ! intent(in)
721 0 : thlp2(i,1) / dt, & ! intent(in)
722 0 : stats_zm(i) ) ! intent(inout)
723 :
724 : call stat_end_update_pt( stats_metadata%irtp2_sf, 1, & ! intent(in)
725 : rtp2(i,1) / dt, & ! intent(in)
726 0 : stats_zm(i) ) ! intent(inout)
727 :
728 : call stat_end_update_pt( stats_metadata%irtpthlp_sf, 1, & ! intent(in)
729 : rtpthlp(i,1) / dt, & ! intent(in)
730 0 : stats_zm(i) ) ! intent(inout)
731 :
732 : call stat_end_update_pt( stats_metadata%iup2_sf, 1, & ! intent(in)
733 : up2(i,1) / dt, & ! intent(in)
734 0 : stats_zm(i) ) ! intent(inout)
735 :
736 : call stat_end_update_pt( stats_metadata%ivp2_sf, 1, & ! intent(in)
737 : vp2(i,1) / dt, & ! intent(in)
738 0 : stats_zm(i) ) ! intent(inout)
739 :
740 : call stat_end_update_pt( stats_metadata%iwp2_sf, 1, & ! intent(in)
741 : wp2(i,1) / dt, & ! intent(in)
742 0 : stats_zm(i) ) ! intent(inout)
743 : end do
744 : end if
745 :
746 352944 : if ( clubb_at_least_debug_level( 2 ) ) then
747 :
748 : !$acc update host( wp2, up2, vp2, thlp2, rtp2, rtpthlp, &
749 : !$acc sclrp2, sclrprtp, sclrpthlp, &
750 : !$acc upwp_sfc, vpwp_sfc, wpthlp, wprtp_sfc )
751 :
752 0 : do i = 1, ngrdcol
753 0 : call sfc_varnce_check( wp2(i,1), up2(i,1), vp2(i,1), & ! intent(in)
754 : thlp2(i,1), rtp2(i,1), rtpthlp(i,1), & ! intent(in)
755 0 : sclrp2(i,1,:), sclrprtp(i,1,:), sclrpthlp(i,1,:) ) ! intent(in)
756 : end do
757 :
758 0 : if ( err_code == clubb_fatal_error ) then
759 :
760 0 : write(fstderr,*) "Error in calc_sfc_varnce"
761 0 : write(fstderr,*) "Intent(in)"
762 :
763 0 : write(fstderr,*) "upwp_sfc = ", upwp_sfc
764 0 : write(fstderr,*) "vpwp_sfc = ", vpwp_sfc
765 0 : write(fstderr,*) "wpthlp = ", wpthlp
766 0 : write(fstderr,*) "wprtp_sfc = ", wprtp_sfc
767 :
768 0 : if ( sclr_dim > 0 ) then
769 0 : write(fstderr,*) "wpsclrp = ", wpsclrp_sfc
770 : endif
771 :
772 0 : write(fstderr,*) "Intent(out)"
773 :
774 0 : write(fstderr,*) "wp2 = ", wp2
775 0 : write(fstderr,*) "up2 = ", up2
776 0 : write(fstderr,*) "vp2 = ", vp2
777 0 : write(fstderr,*) "thlp2 = ", thlp2
778 0 : write(fstderr,*) "rtp2 = ", rtp2
779 0 : write(fstderr,*) "rtpthlp = ", rtpthlp
780 :
781 0 : if ( sclr_dim > 0 ) then
782 0 : write(fstderr,*) "sclrp2 = ", sclrp2
783 0 : write(fstderr,*) "sclrprtp = ", sclrprtp
784 0 : write(fstderr,*) "sclrpthlp = ", sclrpthlp
785 : endif
786 :
787 : endif ! err_code == clubb_fatal_error
788 :
789 : endif ! clubb_at_least_debug_level ( 2 )! Update surface stats
790 :
791 : !$acc exit data delete( uf, depth_pos_wpthlp, min_wp2_sfc_val, &
792 : !$acc um_sfc_sqd, vm_sfc_sqd, usp2_sfc, vsp2_sfc, &
793 : !$acc ustar, zeta, wp2_splat_sfc_correction )
794 :
795 352944 : return
796 :
797 : end subroutine calc_sfc_varnce
798 :
799 : !===============================================================================
800 :
801 : end module sfc_varnce_module
|