Line data Source code
1 :
2 : module advance_windm_edsclrm_module
3 :
4 : implicit none
5 :
6 : private ! Set Default Scope
7 :
8 : public :: advance_windm_edsclrm
9 :
10 : private :: windm_edsclrm_solve, &
11 : compute_uv_tndcy, &
12 : windm_edsclrm_lhs, &
13 : windm_edsclrm_rhs
14 :
15 :
16 : ! Private named constants to avoid string comparisons
17 : integer, parameter, private :: &
18 : windm_edsclrm_um = 1, & ! Named constant to handle um solves
19 : windm_edsclrm_vm = 2, & ! Named constant to handle vm solves
20 : windm_edsclrm_scalar = 3, & ! Named constant to handle scalar solves
21 : clip_upwp = 10, & ! Named constant for upwp clipping
22 : ! NOTE: This must be the same as the clip_upwp
23 : ! declared in clip_explicit!
24 : clip_vpwp = 11 ! Named constant for vpwp clipping
25 : ! NOTE: This must be the same as the clip_vpwp
26 : ! declared in clip_explicit!
27 :
28 : contains
29 :
30 : !=============================================================================
31 352944 : subroutine advance_windm_edsclrm( nz, ngrdcol, gr, dt, &
32 352944 : wm_zt, Km_zm, Kmh_zm, &
33 352944 : ug, vg, um_ref, vm_ref, &
34 352944 : wp2, up2, vp2, um_forcing, vm_forcing, &
35 352944 : edsclrm_forcing, &
36 352944 : rho_ds_zm, invrs_rho_ds_zt, &
37 352944 : fcor, l_implemented, &
38 : nu_vert_res_dep, &
39 : tridiag_solve_method, &
40 : l_predict_upwp_vpwp, &
41 : l_upwind_xm_ma, &
42 : l_uv_nudge, &
43 : l_tke_aniso, &
44 : l_lmm_stepping, &
45 : l_linearize_pbl_winds, &
46 : order_xp2_xpyp, order_wp2_wp3, order_windm, &
47 : stats_metadata, &
48 352944 : stats_zt, stats_zm, stats_sfc, &
49 352944 : um, vm, edsclrm, &
50 352944 : upwp, vpwp, wpedsclrp, &
51 352944 : um_pert, vm_pert, upwp_pert, vpwp_pert )
52 :
53 : ! Description:
54 : ! Solves for both mean horizontal wind components, um and vm, and for the
55 : ! eddy-scalars (passive scalars that don't use the high-order closure).
56 :
57 : ! Uses the LAPACK tridiagonal solver subroutine with 2 + # of scalar(s)
58 : ! back substitutions (since the left hand side matrix is the same for all
59 : ! input variables).
60 :
61 : ! References:
62 : ! Eqn. 8 & 9 on p. 3545 of
63 : ! ``A PDF-Based Model for Boundary Layer Clouds. Part I:
64 : ! Method and Model Description'' Golaz, et al. (2002)
65 : ! JAS, Vol. 59, pp. 3540--3551.
66 : !-----------------------------------------------------------------------
67 :
68 : use grid_class, only: &
69 : grid, & ! Type
70 : zm2zt
71 :
72 : use parameters_model, only: &
73 : ts_nudge, & ! Variable(s)
74 : edsclr_dim
75 :
76 : use parameters_tunable, only: &
77 : nu_vertical_res_dep ! Type(s)
78 :
79 : use clubb_precision, only: &
80 : core_rknd ! Variable(s)
81 :
82 : use stats_type_utilities, only: &
83 : stat_begin_update, & ! Subroutines
84 : stat_end_update, &
85 : stat_update_var
86 :
87 : use stats_variables, only: &
88 : stats_metadata_type
89 :
90 : use clip_explicit, only: &
91 : clip_covar ! Procedure(s)
92 :
93 : use error_code, only: &
94 : clubb_at_least_debug_level, & ! Procedure
95 : err_code, & ! Error Indicator
96 : clubb_fatal_error ! Constant
97 :
98 : use constants_clubb, only: &
99 : one_half, & ! Constant(s)
100 : zero, &
101 : fstderr, &
102 : eps
103 :
104 : use sponge_layer_damping, only: &
105 : uv_sponge_damp_settings, &
106 : uv_sponge_damp_profile, &
107 : sponge_damp_xm ! Procedure(s)
108 :
109 : use stats_type, only: stats ! Type
110 :
111 : use diffusion, only: &
112 : diffusion_zt_lhs ! Procedure(s)
113 :
114 : use mean_adv, only: &
115 : term_ma_zt_lhs ! Procedures
116 :
117 : use model_flags, only: &
118 : l_upwind_Kh_dp_term
119 :
120 : use advance_helper_module, only: &
121 : calc_xpwp
122 :
123 : implicit none
124 :
125 : ! ------------------------ Input Variables ------------------------
126 : integer, intent(in) :: &
127 : nz, &
128 : ngrdcol
129 :
130 : type (grid), target, intent(in) :: &
131 : gr
132 :
133 : real( kind = core_rknd ), intent(in) :: &
134 : dt ! Model timestep [s]
135 :
136 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
137 : wm_zt, & ! w wind component on thermodynamic levels [m/s]
138 : Km_zm, & ! Eddy diffusivity of winds on momentum levs. [m^2/s]
139 : Kmh_zm, & ! Eddy diffusivity of themo on momentum levs. [m^s/s]
140 : ug, & ! u (west-to-east) geostrophic wind comp. [m/s]
141 : vg, & ! v (south-to-north) geostrophic wind comp. [m/s]
142 : um_ref, & ! Reference u wind component for nudging [m/s]
143 : vm_ref, & ! Reference v wind component for nudging [m/s]
144 : wp2, & ! w'^2 (momentum levels) [m^2/s^2]
145 : up2, & ! u'^2 (momentum levels) [m^2/s^2]
146 : vp2, & ! v'^2 (momentum levels) [m^2/s^2]
147 : um_forcing, & ! u forcing [m/s/s]
148 : vm_forcing, & ! v forcing [m/s/s]
149 : rho_ds_zm, & ! Dry, static density on momentum levels [kg/m^3]
150 : invrs_rho_ds_zt ! Inv. dry, static density at thermo. levels [m^3/kg]
151 :
152 : real( kind = core_rknd ), dimension(ngrdcol,nz,edsclr_dim), intent(in) :: &
153 : edsclrm_forcing ! Eddy scalar large-scale forcing [{units vary}/s]
154 :
155 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
156 : fcor ! Coriolis parameter [s^-1]
157 :
158 : logical, intent(in) :: &
159 : l_implemented ! Flag for CLUBB being implemented in a larger model.
160 :
161 : type(nu_vertical_res_dep), intent(in) :: &
162 : nu_vert_res_dep ! Vertical resolution dependent nu values
163 :
164 : integer, intent(in) :: &
165 : tridiag_solve_method ! Specifier for method to solve tridiagonal systems
166 :
167 : logical, intent(in) :: &
168 : l_predict_upwp_vpwp, & ! Flag to predict <u'w'> and <v'w'> along with <u> and <v> alongside
169 : ! the advancement of <rt>, <w'rt'>, <thl>, <wpthlp>, <sclr>, and
170 : ! <w'sclr'> in subroutine advance_xm_wpxp. Otherwise, <u'w'> and
171 : ! <v'w'> are still approximated by eddy diffusivity when <u> and <v>
172 : ! are advanced in subroutine advance_windm_edsclrm.
173 : l_upwind_xm_ma, & ! This flag determines whether we want to use an upwind differencing
174 : ! approximation rather than a centered differencing for turbulent or
175 : ! mean advection terms. It affects rtm, thlm, sclrm, um and vm.
176 : l_uv_nudge, & ! For wind speed nudging
177 : l_tke_aniso, & ! For anisotropic turbulent kinetic energy, i.e. TKE = 1/2
178 : ! (u'^2 + v'^2 + w'^2)
179 : l_lmm_stepping, & ! Apply Linear Multistep Method (LMM) Stepping
180 : l_linearize_pbl_winds ! Flag (used by E3SM) to linearize PBL winds
181 :
182 : integer, intent(in) :: &
183 : order_xp2_xpyp, &
184 : order_wp2_wp3, &
185 : order_windm
186 :
187 : type (stats_metadata_type), intent(in) :: &
188 : stats_metadata
189 :
190 : ! ------------------------ Input/Output Variables ------------------------
191 : type (stats), dimension(ngrdcol), target, intent(inout) :: &
192 : stats_zt, &
193 : stats_zm, &
194 : stats_sfc
195 :
196 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(inout) :: &
197 : um, & ! Mean u (west-to-east) wind component [m/s]
198 : vm, & ! Mean v (south-to-north) wind component [m/s]
199 : upwp, & ! <u'w'> (momentum levels) [m^2/s^2]
200 : vpwp ! <v'w'> (momentum levels) [m^2/s^2]
201 :
202 : ! Input/Output Variable for eddy-scalars
203 : real( kind = core_rknd ), dimension(ngrdcol,nz,edsclr_dim), intent(inout) :: &
204 : edsclrm ! Mean eddy scalar quantity [units vary]
205 :
206 : ! Output Variable for eddy-scalars
207 : real( kind = core_rknd ), dimension(ngrdcol,nz,edsclr_dim), intent(inout) :: &
208 : wpedsclrp ! w'edsclr' (momentum levels) [m/s {units vary}]
209 :
210 : ! Variables used to track perturbed version of winds.
211 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(inout) :: &
212 : um_pert, & ! perturbed <u> [m/s]
213 : vm_pert, & ! perturbed <v> [m/s]
214 : upwp_pert, & ! perturbed <u'w'> [m^2/s^2]
215 : vpwp_pert ! perturbed <v'w'> [m^2/s^2]
216 :
217 : ! ------------------------ Local Variables ------------------------
218 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
219 705888 : um_old, & ! Saved value of mean u (west-to-east) wind component [m/s]
220 705888 : vm_old ! Saved value of Mean v (south-to-north) wind component [m/s]
221 :
222 : real( kind = core_rknd ), dimension(ngrdcol,nz,edsclr_dim) :: &
223 705888 : edsclrm_old ! Saved value of mean eddy scalar quantity [units vary]
224 :
225 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
226 705888 : um_tndcy, & ! u wind component tendency [m/s^2]
227 705888 : vm_tndcy ! v wind component tendency [m/s^2]
228 :
229 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
230 705888 : upwp_chnge, & ! Net change of u'w' due to clipping [m^2/s^2]
231 705888 : vpwp_chnge ! Net change of v'w' due to clipping [m^2/s^2]
232 :
233 : real( kind = core_rknd ), dimension(3,ngrdcol,nz) :: &
234 705888 : lhs ! The implicit part of the tridiagonal matrix [units vary]
235 :
236 : real( kind = core_rknd ), dimension(ngrdcol,nz,max(2,edsclr_dim)) :: &
237 705888 : rhs, &! The explicit part of the tridiagonal matrix [units vary]
238 705888 : solution ! The solution to the tridiagonal matrix [units vary]
239 :
240 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
241 705888 : wind_speed, & ! wind speed; sqrt(u^2 + v^2) [m/s]
242 705888 : wind_speed_pert ! perturbed wind speed; sprt(u^2 + v^2) [m/s]
243 :
244 : real( kind = core_rknd ), dimension(ngrdcol) :: &
245 705888 : u_star_sqd, & ! Surface friction velocity, u_star, squared [m/s]
246 705888 : u_star_sqd_pert ! perturbed u_star, squared [m/s]
247 :
248 : logical :: &
249 : l_imp_sfc_momentum_flux ! Flag for implicit momentum surface fluxes.
250 :
251 : integer :: nrhs ! Number of right hand side terms
252 :
253 : integer :: i, k, edsclr, j ! Array index
254 :
255 : logical :: l_first_clip_ts, l_last_clip_ts ! flags for clip_covar
256 :
257 : real( kind = core_rknd ), dimension(ngrdcol) :: &
258 705888 : nu_zero
259 :
260 : real( kind = core_rknd ), dimension(3,ngrdcol,nz) :: &
261 705888 : lhs_diff, & ! LHS diffustion terms
262 705888 : lhs_ma_zt ! LHS mean advection terms
263 :
264 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
265 705888 : Km_zt, & ! Eddy diffusivity of winds on momentum levs. [m^2/s]
266 705888 : Kmh_zt, & ! Eddy diffusivity of themo on momentum levs. [m^s/s]
267 705888 : Km_zm_p_nu10, & ! Km_zm plus nu_vert_res_dep%nu10
268 705888 : xpwp ! x'w' for arbitrary x
269 :
270 : integer, parameter :: &
271 : kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
272 : k_tdiag = 2, & ! Thermodynamic main diagonal index.
273 : km1_tdiag = 3 ! Thermodynamic subdiagonal index.
274 :
275 : ! Whether perturbed winds are being solved.
276 : logical :: l_perturbed_wind
277 :
278 : ! ------------------------ Begin Code ------------------------
279 :
280 : !$acc enter data create( um_old, vm_old, um_tndcy, vm_tndcy, &
281 : !$acc upwp_chnge, vpwp_chnge, lhs, rhs, solution, wind_speed, &
282 : !$acc wind_speed_pert, u_star_sqd, u_star_sqd_pert, &
283 : !$acc nu_zero, lhs_diff, lhs_ma_zt, Km_zt, Kmh_zt, &
284 : !$acc Km_zm_p_nu10, xpwp )
285 :
286 : !$acc enter data if( edsclr_dim > 0) create( edsclrm_old )
287 :
288 : !$acc parallel loop gang vector default(present)
289 5893344 : do i = 1, ngrdcol
290 5893344 : nu_zero(i) = zero
291 : end do
292 : !$acc end parallel loop
293 :
294 : !$acc parallel loop gang vector collapse(2) default(present)
295 30353184 : do k = 1, nz
296 501287184 : do i = 1, ngrdcol
297 500934240 : Km_zm_p_nu10(i,k) = Km_zm(i,k) + nu_vert_res_dep%nu10(i)
298 : end do
299 : end do
300 : !$acc end parallel loop
301 :
302 352944 : l_perturbed_wind = ( .not. l_predict_upwp_vpwp ) .and. l_linearize_pbl_winds
303 :
304 352944 : if ( .not. l_implemented ) then
305 : call term_ma_zt_lhs( nz, ngrdcol, wm_zt, gr%weights_zt2zm, & ! intent(in)
306 : gr%invrs_dzt, gr%invrs_dzm, & ! intent(in)
307 : l_upwind_xm_ma, & ! intent(in)
308 0 : lhs_ma_zt ) ! intent(out)
309 : else
310 : !$acc parallel loop gang vector collapse(3) default(present)
311 30353184 : do k = 1, nz
312 501287184 : do i = 1, ngrdcol
313 1913736240 : do j = 1, 3
314 1883736000 : lhs_ma_zt(j,i,k) = zero
315 : end do
316 : end do
317 : end do
318 : !$acc end parallel loop
319 : end if
320 :
321 352944 : if ( .not. l_predict_upwp_vpwp ) then
322 :
323 0 : Km_zt(:,:) = zm2zt( nz, ngrdcol, gr, Km_zm(:,:) )
324 :
325 : !$acc parallel loop gang vector collapse(2) default(present)
326 0 : do k = 1, nz
327 0 : do i = 1, ngrdcol
328 0 : Km_zt(i,k) = max( Km_zt(i,k), zero )
329 : end do
330 : end do
331 : !$acc end parallel loop
332 :
333 : ! Calculate diffusion terms
334 : call diffusion_zt_lhs( nz, ngrdcol, gr, Km_zm, Km_zt, nu_vert_res_dep%nu10, & ! In
335 : invrs_rho_ds_zt, rho_ds_zm, & ! In
336 0 : lhs_diff ) ! Out
337 :
338 : ! The lower boundary condition needs to be applied here at level 2.
339 : if ( .not. l_upwind_Kh_dp_term ) then
340 : !$acc parallel loop gang vector default(present)
341 0 : do i = 1, ngrdcol
342 : ! Thermodynamic superdiagonal: [ x xm(k+1,<t+1>) ]
343 0 : lhs_diff(kp1_tdiag,i,2) = - gr%invrs_dzt(i,2) * invrs_rho_ds_zt(i,2) &
344 0 : * ( Km_zm(i,2) + nu_vert_res_dep%nu10(i) ) &
345 0 : * rho_ds_zm(i,2) * gr%invrs_dzm(i,2)
346 :
347 : ! Thermodynamic main diagonal: [ x xm(k,<t+1>) ]
348 0 : lhs_diff(k_tdiag,i,2) = + gr%invrs_dzt(i,2) * invrs_rho_ds_zt(i,2) &
349 : * ( Km_zm(i,2) + nu_vert_res_dep%nu10(i) ) &
350 0 : * rho_ds_zm(i,2) * gr%invrs_dzm(i,2)
351 :
352 : ! Thermodynamic subdiagonal: [ x xm(k-1,<t+1>) ]
353 0 : lhs_diff(km1_tdiag,i,2) = zero
354 : end do
355 : !$acc end parallel loop
356 : else
357 : !$acc parallel loop gang vector default(present)
358 : do i = 1, ngrdcol
359 : ! Thermodynamic superdiagonal: [ x xm(k+1,<t+1>) ]
360 : lhs_diff(kp1_tdiag,i,2) = - gr%invrs_dzt(i,2) &
361 : * ( Km_zt(i,2) + nu_vert_res_dep%nu10(i) ) &
362 : * gr%invrs_dzm(i,2)
363 :
364 : ! Thermodynamic main diagonal: [ x xm(k,<t+1>) ]
365 : lhs_diff(k_tdiag,i,2) = + gr%invrs_dzt(i,2) &
366 : * ( Km_zt(i,2) + nu_vert_res_dep%nu10(i) ) &
367 : * gr%invrs_dzm(i,2)
368 :
369 : ! Thermodynamic subdiagonal: [ x xm(k-1,<t+1>) ]
370 : lhs_diff(km1_tdiag,i,2) = zero
371 : end do
372 : !$acc end parallel loop
373 : end if
374 :
375 0 : if ( l_lmm_stepping ) then
376 : !$acc parallel loop gang vector collapse(2) default(present)
377 0 : do k = 1, nz
378 0 : do i = 1, ngrdcol
379 0 : um_old(i,k) = um(i,k)
380 0 : vm_old(i,k) = vm(i,k)
381 : end do
382 : end do
383 : !$acc end parallel loop
384 : end if ! l_lmm_stepping
385 :
386 : !----------------------------------------------------------------
387 : ! Prepare tridiagonal system for horizontal winds, um and vm
388 : !----------------------------------------------------------------
389 :
390 : ! Compute Coriolis, geostrophic, and other prescribed wind forcings for um.
391 : call compute_uv_tndcy( nz, ngrdcol, windm_edsclrm_um, & ! intent(in)
392 : fcor, vm, vg, & ! intent(in)
393 : um_forcing, l_implemented, & ! intent(in)
394 : stats_metadata, & ! intent(in)
395 : stats_zt, & ! intent(inout)
396 0 : um_tndcy ) ! intent(out)
397 :
398 : ! Compute Coriolis, geostrophic, and other prescribed wind forcings for vm.
399 : call compute_uv_tndcy( nz, ngrdcol, windm_edsclrm_vm, & ! intent(in)
400 : fcor, um, ug, & ! intent(in)
401 : vm_forcing, l_implemented, & ! intent(in)
402 : stats_metadata, & ! intent(in)
403 : stats_zt, & ! intent(inout)
404 0 : vm_tndcy ) ! intent(out)
405 :
406 : ! Momentum surface fluxes, u'w'|_sfc and v'w'|_sfc, are applied through
407 : ! an implicit method, such that:
408 : ! x'w'|_sfc = - ( u_star(t)^2 / wind_speed(t) ) * xm(t+1).
409 0 : l_imp_sfc_momentum_flux = .true.
410 :
411 : ! Compute wind speed (use threshold "eps" to prevent divide-by-zero
412 : ! error).
413 : !$acc parallel loop gang vector collapse(2) default(present)
414 0 : do k = 1, nz
415 0 : do i = 1, ngrdcol
416 0 : wind_speed(i,k) = max( sqrt( um(i,k)**2 + vm(i,k)**2 ), eps )
417 : end do
418 : end do
419 : !$acc end parallel loop
420 :
421 : ! Compute u_star_sqd according to the definition of u_star.
422 : !$acc parallel loop gang vector default(present)
423 0 : do i = 1, ngrdcol
424 0 : u_star_sqd(i) = sqrt( upwp(i,1)**2 + vpwp(i,1)**2 )
425 : end do
426 : !$acc end parallel loop
427 :
428 : ! Compute the explicit portion of the um equation.
429 : ! Build the right-hand side vector.
430 : call windm_edsclrm_rhs( nz, ngrdcol, gr, windm_edsclrm_um, dt, & ! intent(in)
431 : lhs_diff, um, um_tndcy, & ! intent(in)
432 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
433 : l_imp_sfc_momentum_flux, upwp(:,1), & ! intent(in)
434 : stats_metadata, & ! intent(in)
435 : stats_zt, & ! intent(inout)
436 0 : rhs(:,:,windm_edsclrm_um) ) ! intent(out)
437 :
438 : ! Compute the explicit portion of the vm equation.
439 : ! Build the right-hand side vector.
440 : call windm_edsclrm_rhs( nz, ngrdcol, gr, windm_edsclrm_vm, dt, & ! intent(in)
441 : lhs_diff, vm, vm_tndcy, & ! intent(in)
442 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
443 : l_imp_sfc_momentum_flux, vpwp(:,1), & ! intent(in)
444 : stats_metadata, & ! intent(in)
445 : stats_zt, & ! intent(inout)
446 0 : rhs(:,:,windm_edsclrm_vm) ) ! intent(out)
447 :
448 : ! Store momentum flux (explicit component)
449 :
450 : ! The surface flux, x'w'(1) = x'w'|_sfc, is set elsewhere in the model.
451 : ! upwp(1) = upwp_sfc
452 : ! vpwp(1) = vpwp_sfc
453 :
454 : call calc_xpwp( nz, ngrdcol, gr, &
455 : Km_zm_p_nu10, um, &
456 0 : xpwp )
457 :
458 : ! Solve for x'w' at all intermediate model levels.
459 : ! A Crank-Nicholson timestep is used.
460 : !$acc parallel loop gang vector collapse(2) default(present)
461 0 : do k = 2, nz-1
462 0 : do i = 1, ngrdcol
463 0 : upwp(i,k) = -one_half * xpwp(i,k)
464 : end do
465 : end do
466 : !$acc end parallel loop
467 :
468 : call calc_xpwp( nz, ngrdcol, gr, &
469 : Km_zm_p_nu10, vm, &
470 0 : xpwp )
471 :
472 : !$acc parallel loop gang vector collapse(2) default(present)
473 0 : do k = 2, nz-1
474 0 : do i = 1, ngrdcol
475 0 : vpwp(i,k) = -one_half * xpwp(i,k)
476 : end do
477 : end do
478 : !$acc end parallel loop
479 :
480 : ! A zero-flux boundary condition at the top of the model, d(xm)/dz = 0,
481 : ! means that x'w' at the top model level is 0,
482 : ! since x'w' = - K_zm * d(xm)/dz.
483 : !$acc parallel loop gang vector default(present)
484 0 : do i = 1, ngrdcol
485 0 : upwp(i,nz) = zero
486 0 : vpwp(i,nz) = zero
487 : end do
488 : !$acc end parallel loop
489 :
490 : ! Compute the implicit portion of the um and vm equations.
491 : ! Build the left-hand side matrix.
492 : call windm_edsclrm_lhs( nz, ngrdcol, gr, dt, & ! intent(in)
493 : lhs_ma_zt, lhs_diff, & ! intent(in)
494 : wind_speed, u_star_sqd, & ! intent(in)
495 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
496 : l_implemented, l_imp_sfc_momentum_flux, & ! intent(in)
497 0 : lhs ) ! intent(out)
498 :
499 : ! Decompose and back substitute for um and vm
500 0 : nrhs = 2
501 : call windm_edsclrm_solve( nz, ngrdcol, gr, nrhs, stats_metadata%iwindm_matrix_condt_num, & ! intent(in)
502 : tridiag_solve_method, & ! intent(in)
503 : stats_metadata, & ! intent(in)
504 : stats_sfc, & ! intent(inout)
505 : lhs, rhs, & ! intent(inout)
506 0 : solution ) ! intent(out)
507 :
508 : ! Check for singular matrices and bad LAPACK arguments
509 0 : if ( clubb_at_least_debug_level( 0 ) ) then
510 0 : if ( err_code == clubb_fatal_error ) then
511 0 : write(fstderr,*) "Fatal error solving for um/vm"
512 0 : return
513 : end if
514 : end if
515 :
516 : !----------------------------------------------------------------
517 : ! Update zonal (west-to-east) component of mean wind, um
518 : !----------------------------------------------------------------
519 : !$acc parallel loop gang vector collapse(2) default(present)
520 0 : do k = 1, nz
521 0 : do i = 1, ngrdcol
522 0 : um(i,k) = solution(i,k,windm_edsclrm_um)
523 : end do
524 : end do
525 : !$acc end parallel loop
526 :
527 : !----------------------------------------------------------------
528 : ! Update meridional (south-to-north) component of mean wind, vm
529 : !----------------------------------------------------------------
530 : !$acc parallel loop gang vector collapse(2) default(present)
531 0 : do k = 1, nz
532 0 : do i = 1, ngrdcol
533 0 : vm(i,k) = solution(i,k,windm_edsclrm_vm)
534 : end do
535 : end do
536 : !$acc end parallel loop
537 :
538 0 : if ( stats_metadata%l_stats_samp ) then
539 :
540 : !$acc update host( um, lhs_diff, lhs_ma_zt, invrs_rho_ds_zt, u_star_sqd, &
541 : !$acc rho_ds_zm, wind_speed, vm )
542 :
543 0 : do i = 1, ngrdcol
544 : ! Implicit contributions to um and vm
545 0 : call windm_edsclrm_implicit_stats( nz, windm_edsclrm_um, um(i,:), & ! intent(in)
546 0 : gr%invrs_dzt(i,:), & ! intent(in)
547 0 : lhs_diff(:,i,:), lhs_ma_zt(:,i,:), & ! intent(in)
548 0 : invrs_rho_ds_zt(i,:), u_star_sqd(i), & ! intent(in)
549 : rho_ds_zm(i,:), wind_speed(i,:), & ! intent(in)
550 : l_imp_sfc_momentum_flux, & ! intent(in)
551 : stats_metadata, & ! intent(in)
552 0 : stats_zt(i) ) ! intent(inout)
553 :
554 : call windm_edsclrm_implicit_stats( nz, windm_edsclrm_vm, vm(i,:), & ! intent(in)
555 0 : gr%invrs_dzt(i,:), & ! intent(in)
556 : lhs_diff(:,i,:), lhs_ma_zt(:,i,:), & ! intent(in)
557 : invrs_rho_ds_zt(i,:), u_star_sqd(i), & ! intent(in)
558 : rho_ds_zm(i,:), wind_speed(i,:), & ! intent(in)
559 : l_imp_sfc_momentum_flux, & ! intent(in)
560 : stats_metadata, & ! intent(in)
561 0 : stats_zt(i) ) ! intent(inout)
562 : end do
563 : end if ! stats_metadata%l_stats_samp
564 :
565 : ! The values of um(1) and vm(1) are located below the model surface and
566 : ! do not affect the rest of the model. The values of um(1) or vm(1) are
567 : ! simply set to the values of um(2) and vm(2), respectively, after the
568 : ! equation matrices has been solved. Even though um and vm would sharply
569 : ! decrease to a value of 0 at the surface, this is done to avoid
570 : ! confusion on plots of the vertical profiles of um and vm.
571 : !$acc parallel loop gang vector default(present)
572 0 : do i = 1, ngrdcol
573 0 : um(i,1) = um(i,2)
574 0 : vm(i,1) = vm(i,2)
575 : end do
576 : !$acc end parallel loop
577 :
578 0 : if ( l_lmm_stepping ) then
579 : !$acc parallel loop gang vector collapse(2) default(present)
580 0 : do k = 1, nz
581 0 : do i = 1, ngrdcol
582 0 : um(i,k) = one_half * ( um_old(i,k) + um(i,k) )
583 0 : vm(i,k) = one_half * ( vm_old(i,k) + vm(i,k) )
584 : end do
585 : end do
586 : !$acc end parallel loop
587 : endif ! l_lmm_stepping ) then
588 :
589 0 : if ( uv_sponge_damp_settings%l_sponge_damping ) then
590 :
591 : !$acc update host( um, vm, um_ref, vm_ref )
592 :
593 : ! _sponge_damp_settings and _sponge_damp_profile could potentially vary
594 : ! from column to column, but there is no column index in those variables.
595 : ! Thus this code is potentially unsafe when implemented in a host model,
596 : ! which is indicated by l_implemented = T
597 0 : if ( l_implemented ) then
598 : write(fstderr,*) "l_sponge_damping = T and l_implemented = T &
599 0 : -- this is likely unsafe and considered fatal"
600 0 : err_code = clubb_fatal_error
601 0 : return
602 : end if
603 :
604 0 : if ( stats_metadata%l_stats_samp ) then
605 0 : do i = 1, ngrdcol
606 0 : call stat_begin_update( nz, stats_metadata%ium_sdmp, um(i,:) / dt, & ! intent(in)
607 0 : stats_zt(i) ) ! intent(inout)
608 : call stat_begin_update( nz, stats_metadata%ivm_sdmp, vm(i,:) / dt, & ! intent(in)
609 0 : stats_zt(i) ) ! intent(inout)
610 : end do
611 : end if
612 :
613 0 : do i = 1, ngrdcol
614 0 : um(i,1:nz) = sponge_damp_xm( nz, dt, gr%zt(i,:), gr%zm(i,:), &
615 0 : um_ref(i,1:nz), um(i,1:nz), uv_sponge_damp_profile )
616 : end do
617 :
618 0 : do i = 1, ngrdcol
619 0 : vm(i,1:nz) = sponge_damp_xm( nz, dt, gr%zt(i,:), gr%zm(i,:), &
620 0 : vm_ref(i,1:nz), vm(i,1:nz), uv_sponge_damp_profile )
621 : end do
622 :
623 0 : if ( stats_metadata%l_stats_samp ) then
624 0 : do i = 1, ngrdcol
625 0 : call stat_end_update( nz, stats_metadata%ium_sdmp, um(i,:) / dt, & ! intent(in)
626 0 : stats_zt(i) ) ! intent(inout)
627 : call stat_end_update( nz, stats_metadata%ivm_sdmp, vm(i,:) / dt, & ! intent(in)
628 0 : stats_zt(i) ) ! intent(inout)
629 : end do
630 : end if
631 :
632 : !$acc update device( um, vm )
633 :
634 : end if ! uv_sponge_damp_settings%l_sponge_damping
635 :
636 : ! Second part of momentum (implicit component)
637 :
638 : ! Solve for x'w' at all intermediate model levels.
639 : ! A Crank-Nicholson timestep is used.
640 : call calc_xpwp( nz, ngrdcol, gr, &
641 : Km_zm_p_nu10, um, &
642 0 : xpwp )
643 :
644 : !$acc parallel loop gang vector collapse(2) default(present)
645 0 : do k = 2, nz-1
646 0 : do i = 1, ngrdcol
647 0 : upwp(i,k) = upwp(i,k) - one_half * xpwp(i,k)
648 : end do
649 : end do
650 : !$acc end parallel loop
651 :
652 : call calc_xpwp( nz, ngrdcol, gr, &
653 : Km_zm_p_nu10, vm, &
654 0 : xpwp )
655 :
656 : !$acc parallel loop gang vector collapse(2) default(present)
657 0 : do k = 2, nz-1
658 0 : do i = 1, ngrdcol
659 0 : vpwp(i,k) = vpwp(i,k) - one_half * xpwp(i,k)
660 : end do
661 : end do
662 : !$acc end parallel loop
663 :
664 : ! Adjust um and vm if nudging is turned on.
665 0 : if ( l_uv_nudge ) then
666 :
667 : ! Reflect nudging in budget
668 0 : if ( stats_metadata%l_stats_samp ) then
669 : !$acc update host( um, vm )
670 0 : do i = 1, ngrdcol
671 0 : call stat_begin_update( nz, stats_metadata%ium_ndg, um(i,:) / dt, & ! intent(in)
672 0 : stats_zt(i) ) ! intent(inout)
673 : call stat_begin_update( nz, stats_metadata%ivm_ndg, vm(i,:) / dt, & ! intent(in)
674 0 : stats_zt(i) ) ! intent(inout)
675 : end do
676 : end if
677 :
678 : !$acc parallel loop gang vector collapse(2) default(present)
679 0 : do k = 1, nz
680 0 : do i = 1, ngrdcol
681 0 : um(i,k) = um(i,k) - ( ( um(i,k) - um_ref(i,k) ) * (dt/ts_nudge) )
682 0 : vm(i,k) = vm(i,k) - ( ( vm(i,k) - vm_ref(i,k) ) * (dt/ts_nudge) )
683 : end do
684 : end do
685 : !$acc end parallel loop
686 :
687 0 : if ( stats_metadata%l_stats_samp ) then
688 : !$acc update host( um, vm )
689 0 : do i = 1, ngrdcol
690 0 : call stat_end_update( nz, stats_metadata%ium_ndg, um(i,:) / dt, & ! intent(in)
691 0 : stats_zt(i) ) ! intent(inout)
692 : call stat_end_update( nz, stats_metadata%ivm_ndg, vm(i,:) / dt, & ! intent(in)
693 0 : stats_zt(i) ) ! intent(inout)
694 : end do
695 : end if
696 :
697 : end if ! l_uv_nudge
698 :
699 0 : if ( stats_metadata%l_stats_samp ) then
700 : !$acc update host( um_ref, vm_ref )
701 0 : do i = 1, ngrdcol
702 0 : call stat_update_var( stats_metadata%ium_ref, um_ref(i,:), & ! intent(in)
703 0 : stats_zt(i) ) ! intent(inout)
704 : call stat_update_var( stats_metadata%ivm_ref, vm_ref(i,:), & ! intent(in)
705 0 : stats_zt(i) ) ! intent(inout)
706 : end do
707 : end if
708 :
709 : if ( order_windm < order_wp2_wp3 &
710 0 : .and. order_windm < order_xp2_xpyp ) then
711 0 : l_first_clip_ts = .true.
712 0 : l_last_clip_ts = .false.
713 : elseif ( order_windm > order_wp2_wp3 &
714 0 : .and. order_windm > order_xp2_xpyp ) then
715 0 : l_first_clip_ts = .false.
716 0 : l_last_clip_ts = .true.
717 : else
718 0 : l_first_clip_ts = .false.
719 0 : l_last_clip_ts = .false.
720 : endif
721 :
722 0 : if ( l_tke_aniso ) then
723 :
724 : ! Clipping for u'w'
725 : !
726 : ! Clipping u'w' at each vertical level, based on the
727 : ! correlation of u and w at each vertical level, such that:
728 : ! corr_(u,w) = u'w' / [ sqrt(u'^2) * sqrt(w'^2) ];
729 : ! -1 <= corr_(u,w) <= 1.
730 : !
731 : ! Since u'^2, w'^2, and u'w' are each advanced in different
732 : ! subroutines from each other in advance_clubb_core, clipping for u'w'
733 : ! has to be done three times during each timestep (once after each
734 : ! variable has been updated).
735 : ! This is the third instance of u'w' clipping.
736 : !l_first_clip_ts = .false.
737 : !l_last_clip_ts = .true.
738 : call clip_covar( nz, ngrdcol, gr, clip_upwp, l_first_clip_ts, & ! intent(in)
739 : l_last_clip_ts, dt, wp2, up2, & ! intent(in)
740 : l_predict_upwp_vpwp, & ! intent(in)
741 : stats_metadata, & ! intent(in)
742 : stats_zm, & ! intent(inout)
743 0 : upwp, upwp_chnge ) ! intent(inout)
744 :
745 : ! Clipping for v'w'
746 : !
747 : ! Clipping v'w' at each vertical level, based on the
748 : ! correlation of v and w at each vertical level, such that:
749 : ! corr_(v,w) = v'w' / [ sqrt(v'^2) * sqrt(w'^2) ];
750 : ! -1 <= corr_(v,w) <= 1.
751 : !
752 : ! Since v'^2, w'^2, and v'w' are each advanced in different
753 : ! subroutines from each other in advance_clubb_core, clipping for v'w'
754 : ! has to be done three times during each timestep (once after each
755 : ! variable has been updated).
756 : ! This is the third instance of v'w' clipping.
757 : !l_first_clip_ts = .false.
758 : !l_last_clip_ts = .true.
759 : call clip_covar( nz, ngrdcol, gr, clip_vpwp, l_first_clip_ts, & ! intent(in)
760 : l_last_clip_ts, dt, wp2, vp2, & ! intent(in)
761 : l_predict_upwp_vpwp, & ! intent(in)
762 : stats_metadata, & ! intent(in)
763 : stats_zm, & ! intent(inout)
764 0 : vpwp, vpwp_chnge ) ! intent(inout)
765 : else
766 :
767 : ! intent(in) this case, it is assumed that
768 : ! u'^2 == v'^2 == w'^2, and the variables `up2' and `vp2' do not
769 : ! interact with any other variables.
770 : !l_first_clip_ts = .false.
771 : !l_last_clip_ts = .true.
772 : call clip_covar( nz, ngrdcol, gr, clip_upwp, l_first_clip_ts, & ! intent(in)
773 : l_last_clip_ts, dt, wp2, wp2, & ! intent(in)
774 : l_predict_upwp_vpwp, & ! intent(in)
775 : stats_metadata, & ! intent(in)
776 : stats_zm, & ! intent(inout)
777 0 : upwp, upwp_chnge ) ! intent(inout)
778 :
779 : call clip_covar( nz, ngrdcol, gr, clip_vpwp, l_first_clip_ts, & ! intent(in)
780 : l_last_clip_ts, dt, wp2, wp2, & ! intent(in)
781 : l_predict_upwp_vpwp, & ! intent(in)
782 : stats_metadata, & ! intent(in)
783 : stats_zm, & ! intent(inout)
784 0 : vpwp, vpwp_chnge ) ! intent(inout)
785 : endif ! l_tke_aniso
786 :
787 : endif ! .not. l_predict_upwp_vpwp
788 :
789 :
790 352944 : if ( l_perturbed_wind ) then
791 :
792 : !----------------------------------------------------------------
793 : ! Prepare tridiagonal system for horizontal winds, um and vm
794 : !----------------------------------------------------------------
795 :
796 : ! Momentum surface fluxes, u'w'|_sfc and v'w'|_sfc, are applied through
797 : ! an implicit method, such that:
798 : ! x'w'|_sfc = - ( u_star(t)^2 / wind_speed(t) ) * xm(t+1).
799 0 : l_imp_sfc_momentum_flux = .true.
800 :
801 : ! Compute wind speed (use threshold "eps" to prevent divide-by-zero
802 : ! error).
803 : !$acc parallel loop gang vector collapse(2) default(present)
804 0 : do k = 1, nz
805 0 : do i = 1, ngrdcol
806 0 : wind_speed_pert(i,k) = max( sqrt( (um_pert(i,k))**2 + (vm_pert(i,k))**2 ), eps )
807 : end do
808 : end do
809 : !$acc end parallel loop
810 :
811 : ! Compute u_star_sqd according to the definition of u_star.
812 : !$acc parallel loop gang vector default(present)
813 0 : do i = 1, ngrdcol
814 0 : u_star_sqd_pert(i) = sqrt( upwp_pert(i,1)**2 + vpwp_pert(i,1)**2 )
815 : end do
816 : !$acc end parallel loop
817 :
818 : ! Compute the explicit portion of the um equation.
819 : ! Build the right-hand side vector.
820 : call windm_edsclrm_rhs( nz, ngrdcol, gr, windm_edsclrm_um, dt, & ! intent(in)
821 : lhs_diff, um_pert, um_tndcy, & ! intent(in)
822 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
823 : l_imp_sfc_momentum_flux, upwp_pert(:,1), & ! intent(in)
824 : stats_metadata, & ! intent(in)
825 : stats_zt, & ! intent(inout)
826 0 : rhs(:,:,windm_edsclrm_um) ) ! intent(out)
827 :
828 : ! Compute the explicit portion of the vm equation.
829 : ! Build the right-hand side vector.
830 : call windm_edsclrm_rhs( nz, ngrdcol, gr, windm_edsclrm_vm, dt, & ! intent(in)
831 : lhs_diff, vm_pert, vm_tndcy, & ! intent(in)
832 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
833 : l_imp_sfc_momentum_flux, vpwp_pert(:,1), & ! intent(in)
834 : stats_metadata, & ! intent(in)
835 : stats_zt, & ! intent(inout)
836 0 : rhs(:,:,windm_edsclrm_vm) ) ! intent(out)
837 :
838 : ! Store momentum flux (explicit component)
839 :
840 : ! The surface flux, x'w'(1) = x'w'|_sfc, is set elsewhere in the model.
841 : ! upwp(1) = upwp_sfc
842 : ! vpwp(1) = vpwp_sfc
843 :
844 : ! Solve for x'w' at all intermediate model levels.
845 : ! A Crank-Nicholson timestep is used.
846 : call calc_xpwp( nz, ngrdcol, gr, &
847 : Km_zm_p_nu10, um_pert, &
848 0 : xpwp )
849 :
850 : !$acc parallel loop gang vector collapse(2) default(present)
851 0 : do k = 2, nz-1
852 0 : do i = 1, ngrdcol
853 0 : upwp_pert(i,k) = -one_half * xpwp(i,k)
854 : end do
855 : end do
856 : !$acc end parallel loop
857 :
858 : call calc_xpwp( nz, ngrdcol, gr, &
859 : Km_zm_p_nu10, vm_pert, &
860 0 : xpwp )
861 :
862 : !$acc parallel loop gang vector collapse(2) default(present)
863 0 : do k = 2, nz-1
864 0 : do i = 1, ngrdcol
865 0 : vpwp_pert(i,k) = -one_half * xpwp(i,k)
866 : end do
867 : end do
868 : !$acc end parallel loop
869 :
870 : ! A zero-flux boundary condition at the top of the model, d(xm)/dz = 0,
871 : ! means that x'w' at the top model level is 0,
872 : ! since x'w' = - K_zm * d(xm)/dz.
873 : !$acc parallel loop gang vector default(present)
874 0 : do i = 1, ngrdcol
875 0 : upwp_pert(i,nz) = zero
876 0 : vpwp_pert(i,nz) = zero
877 : end do
878 : !$acc end parallel loop
879 :
880 : ! Compute the implicit portion of the um and vm equations.
881 : ! Build the left-hand side matrix.
882 : call windm_edsclrm_lhs( nz, ngrdcol, gr, dt, & ! intent(in)
883 : lhs_ma_zt, lhs_diff, & ! intent(in)
884 : wind_speed_pert, u_star_sqd_pert, & ! intent(in)
885 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
886 : l_implemented, l_imp_sfc_momentum_flux, & ! intent(in)
887 0 : lhs ) ! intent(out)
888 :
889 : ! Decompose and back substitute for um and vm
890 0 : nrhs = 2
891 : call windm_edsclrm_solve( nz, ngrdcol, gr, nrhs, stats_metadata%iwindm_matrix_condt_num, & ! intent(in)
892 : tridiag_solve_method, & ! intent(in)
893 : stats_metadata, & ! intent(in)
894 : stats_sfc, & ! intent(in)
895 : lhs, rhs, & ! intent(inout)
896 0 : solution ) ! intent(out)
897 :
898 : ! Check for singular matrices and bad LAPACK arguments
899 0 : if ( clubb_at_least_debug_level( 0 ) ) then
900 0 : if ( err_code == clubb_fatal_error ) then
901 0 : write(fstderr,*) "Fatal error solving for um_pert/vm_pert"
902 0 : return
903 : endif
904 : endif
905 :
906 : !----------------------------------------------------------------
907 : ! Update zonal (west-to-east) component of mean wind, um
908 : !----------------------------------------------------------------
909 : !$acc parallel loop gang vector collapse(2) default(present)
910 0 : do k = 1, nz
911 0 : do i = 1, ngrdcol
912 0 : um_pert(i,k) = solution(i,k,windm_edsclrm_um)
913 : end do
914 : end do
915 : !$acc end parallel loop
916 :
917 : !----------------------------------------------------------------
918 : ! Update meridional (south-to-north) component of mean wind, vm
919 : !----------------------------------------------------------------
920 : !$acc parallel loop gang vector collapse(2) default(present)
921 0 : do k = 1, nz
922 0 : do i = 1, ngrdcol
923 0 : vm_pert(i,k) = solution(i,k,windm_edsclrm_vm)
924 : end do
925 : end do
926 : !$acc end parallel loop
927 :
928 : ! The values of um(1) and vm(1) are located below the model surface and
929 : ! do not affect the rest of the model. The values of um(1) or vm(1) are
930 : ! simply set to the values of um(2) and vm(2), respectively, after the
931 : ! equation matrices has been solved. Even though um and vm would sharply
932 : ! decrease to a value of 0 at the surface, this is done to avoid
933 : ! confusion on plots of the vertical profiles of um and vm.
934 : !$acc parallel loop gang vector default(present)
935 0 : do i = 1, ngrdcol
936 0 : um_pert(i,1) = um_pert(i,2)
937 0 : vm_pert(i,1) = vm_pert(i,2)
938 : end do
939 : !$acc end parallel loop
940 :
941 : ! Second part of momentum (implicit component)
942 :
943 : ! Solve for x'w' at all intermediate model levels.
944 : ! A Crank-Nicholson timestep is used.
945 : call calc_xpwp( nz, ngrdcol, gr, &
946 : Km_zm_p_nu10, um_pert, &
947 0 : xpwp )
948 :
949 : !$acc parallel loop gang vector collapse(2) default(present)
950 0 : do k = 2, nz-1
951 0 : do i = 1, ngrdcol
952 0 : upwp_pert(i,k) = upwp_pert(i,k) - one_half * xpwp(i,k)
953 : end do
954 : end do
955 : !$acc end parallel loop
956 :
957 : call calc_xpwp( nz, ngrdcol, gr, &
958 : Km_zm_p_nu10, vm_pert, &
959 0 : xpwp )
960 :
961 : !$acc parallel loop gang vector collapse(2) default(present)
962 0 : do k = 2, nz-1
963 0 : do i = 1, ngrdcol
964 0 : vpwp_pert(i,k) = vpwp_pert(i,k) - one_half * xpwp(i,k)
965 : end do
966 : end do
967 : !$acc end parallel loop
968 :
969 0 : if ( l_tke_aniso ) then
970 :
971 : ! Clipping for u'w'
972 : !
973 : ! Clipping u'w' at each vertical level, based on the
974 : ! correlation of u and w at each vertical level, such that:
975 : ! corr_(u,w) = u'w' / [ sqrt(u'^2) * sqrt(w'^2) ];
976 : ! -1 <= corr_(u,w) <= 1.
977 : !
978 : ! Since u'^2, w'^2, and u'w' are each advanced in different
979 : ! subroutines from each other in advance_clubb_core, clipping for u'w'
980 : ! has to be done three times during each timestep (once after each
981 : ! variable has been updated).
982 : ! This is the third instance of u'w' clipping.
983 0 : l_first_clip_ts = .false.
984 0 : l_last_clip_ts = .true.
985 : call clip_covar( nz, ngrdcol, gr, clip_upwp, l_first_clip_ts, & ! intent(in)
986 : l_last_clip_ts, dt, wp2, up2, & ! intent(in)
987 : l_predict_upwp_vpwp, & ! intent(in)
988 : stats_metadata, & ! intent(in)
989 : stats_zm, & ! intent(inout)
990 0 : upwp_pert, upwp_chnge ) ! intent(inout)
991 :
992 : ! Clipping for v'w'
993 : !
994 : ! Clipping v'w' at each vertical level, based on the
995 : ! correlation of v and w at each vertical level, such that:
996 : ! corr_(v,w) = v'w' / [ sqrt(v'^2) * sqrt(w'^2) ];
997 : ! -1 <= corr_(v,w) <= 1.
998 : !
999 : ! Since v'^2, w'^2, and v'w' are each advanced in different
1000 : ! subroutines from each other in advance_clubb_core, clipping for v'w'
1001 : ! has to be done three times during each timestep (once after each
1002 : ! variable has been updated).
1003 : ! This is the third instance of v'w' clipping.
1004 0 : l_first_clip_ts = .false.
1005 0 : l_last_clip_ts = .true.
1006 : call clip_covar( nz, ngrdcol, gr, clip_vpwp, l_first_clip_ts, & ! intent(in)
1007 : l_last_clip_ts, dt, wp2, vp2, & ! intent(in)
1008 : l_predict_upwp_vpwp, & ! intent(in)
1009 : stats_metadata, & ! intent(in)
1010 : stats_zm, & ! intent(inout)
1011 0 : vpwp_pert, vpwp_chnge ) ! intent(inout)
1012 : else
1013 :
1014 : ! intent(in) this case, it is assumed that
1015 : ! u'^2 == v'^2 == w'^2, and the variables `up2' and `vp2' do not
1016 : ! interact with any other variables.
1017 0 : l_first_clip_ts = .false.
1018 0 : l_last_clip_ts = .true.
1019 : call clip_covar( nz, ngrdcol, gr, clip_upwp, l_first_clip_ts, & ! intent(in)
1020 : l_last_clip_ts, dt, wp2, wp2, & ! intent(in)
1021 : l_predict_upwp_vpwp, & ! intent(in)
1022 : stats_metadata, & ! intent(in)
1023 : stats_zm, & ! intent(inout)
1024 0 : upwp_pert, upwp_chnge ) ! intent(inout)
1025 :
1026 : call clip_covar( nz, ngrdcol, gr, clip_vpwp, l_first_clip_ts, & ! intent(in)
1027 : l_last_clip_ts, dt, wp2, wp2, & ! intent(in)
1028 : l_predict_upwp_vpwp, & ! intent(in)
1029 : stats_metadata, & ! intent(in)
1030 : stats_zm, & ! intent(inout)
1031 0 : vpwp_pert, vpwp_chnge ) ! intent(inout)
1032 :
1033 : end if ! l_tke_aniso
1034 : end if ! l_perturbed_wind
1035 :
1036 : !----------------------------------------------------------------
1037 : ! Prepare tridiagonal system for eddy-scalars
1038 : !----------------------------------------------------------------
1039 :
1040 352944 : if ( edsclr_dim > 0 ) then
1041 :
1042 352944 : Kmh_zt(:,:) = zm2zt( nz, ngrdcol, gr, Kmh_zm(:,:) )
1043 :
1044 : !$acc parallel loop gang vector collapse(2) default(present)
1045 29647296 : do k = 2, nz-1
1046 489500496 : do i = 1, ngrdcol
1047 489147552 : Kmh_zt(i,k) = max( Kmh_zt(i,k), zero )
1048 : end do
1049 : end do
1050 : !$acc end parallel loop
1051 :
1052 : ! Calculate diffusion terms
1053 : call diffusion_zt_lhs( nz, ngrdcol, gr, Kmh_zm, Kmh_zt, nu_zero, & ! intent(in)
1054 : invrs_rho_ds_zt, rho_ds_zm, & ! intent(in)
1055 352944 : lhs_diff ) ! intent(out)
1056 :
1057 : ! The lower boundary condition needs to be applied here at level 2.
1058 : if ( .not. l_upwind_Kh_dp_term ) then
1059 :
1060 : !$acc parallel loop gang vector default(present)
1061 5893344 : do i = 1, ngrdcol
1062 : ! Thermodynamic superdiagonal: [ x xm(k+1,<t+1>) ]
1063 11080800 : lhs_diff(kp1_tdiag,i,2) = - gr%invrs_dzt(i,2) * invrs_rho_ds_zt(i,2) &
1064 : * ( Kmh_zm(i,2) + nu_zero(i) ) &
1065 11080800 : * rho_ds_zm(i,2) * gr%invrs_dzm(i,2)
1066 :
1067 : ! Thermodynamic main diagonal: [ x xm(k,<t+1>) ]
1068 0 : lhs_diff(k_tdiag,i,2) = + gr%invrs_dzt(i,2) * invrs_rho_ds_zt(i,2) &
1069 : * ( Kmh_zm(i,2) + nu_zero(i) ) &
1070 5540400 : * rho_ds_zm(i,2) * gr%invrs_dzm(i,2)
1071 :
1072 : ! Thermodynamic subdiagonal: [ x xm(k-1,<t+1>) ]
1073 5893344 : lhs_diff(km1_tdiag,i,2) = zero
1074 : end do
1075 : !$acc end parallel loop
1076 :
1077 : else
1078 :
1079 : !$acc parallel loop gang vector default(present)
1080 : do i = 1, ngrdcol
1081 : ! Thermodynamic superdiagonal: [ x xm(k+1,<t+1>) ]
1082 : lhs_diff(kp1_tdiag,i,2) = - gr%invrs_dzt(i,2) &
1083 : * ( Kmh_zt(i,2) + nu_zero(i) ) * gr%invrs_dzm(i,2)
1084 :
1085 : ! Thermodynamic main diagonal: [ x xm(k,<t+1>) ]
1086 : lhs_diff(k_tdiag,i,2) = + gr%invrs_dzt(i,2) &
1087 : * ( Kmh_zt(i,2) + nu_zero(i) ) * gr%invrs_dzm(i,2)
1088 :
1089 : ! Thermodynamic subdiagonal: [ x xm(k-1,<t+1>) ]
1090 : lhs_diff(km1_tdiag,i,2) = zero
1091 : end do
1092 : !$acc end parallel loop
1093 :
1094 : end if
1095 :
1096 352944 : if ( l_lmm_stepping ) then
1097 : !$acc parallel loop gang vector collapse(3) default(present)
1098 0 : do j = 1, edsclr_dim
1099 0 : do k = 1, nz
1100 0 : do i = 1, ngrdcol
1101 0 : edsclrm_old(i,k,j) = edsclrm(i,k,j)
1102 : end do
1103 : end do
1104 : end do
1105 : !$acc end parallel loop
1106 : endif ! l_lmm_stepping
1107 :
1108 : ! Eddy-scalar surface fluxes, x'w'|_sfc, are applied through an explicit
1109 : ! method.
1110 352944 : l_imp_sfc_momentum_flux = .false.
1111 :
1112 : ! Compute the explicit portion of eddy scalar equation.
1113 : ! Build the right-hand side vector.
1114 : ! Because of statistics, we have to use a DO rather than a FORALL here
1115 : ! -dschanen 7 Oct 2008
1116 8470656 : do edsclr = 1, edsclr_dim
1117 : call windm_edsclrm_rhs( nz, ngrdcol, gr, windm_edsclrm_scalar, dt, & ! intent(in)
1118 : lhs_diff, edsclrm(:,:,edsclr), & ! intent(in)
1119 : edsclrm_forcing(:,:,edsclr), & ! intent(in)
1120 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
1121 : l_imp_sfc_momentum_flux, wpedsclrp(:,1,edsclr), & ! intent(in)
1122 : stats_metadata, & ! intent(in)
1123 : stats_zt, & ! intent(inout)
1124 8470656 : rhs(:,:,edsclr) ) ! intent(out)
1125 : enddo
1126 :
1127 :
1128 : ! Store momentum flux (explicit component)
1129 :
1130 : ! The surface flux, x'w'(1) = x'w'|_sfc, is set elsewhere in the model.
1131 : ! wpedsclrp(1,1:edsclr_dim) = wpedsclrp_sfc(1:edsclr_dim)
1132 :
1133 : ! Solve for x'w' at all intermediate model levels.
1134 : ! A Crank-Nicholson timestep is used
1135 8470656 : do edsclr = 1, edsclr_dim
1136 :
1137 : call calc_xpwp( nz, ngrdcol, gr, &
1138 : Km_zm_p_nu10, edsclrm(:,:,edsclr), &
1139 8117712 : xpwp )
1140 :
1141 : !$acc parallel loop gang vector collapse(2) default(present)
1142 682240752 : do k = 2, nz-1
1143 11258511408 : do i = 1, ngrdcol
1144 11250393696 : wpedsclrp(i,k,edsclr) = -one_half * xpwp(i,k)
1145 : end do
1146 : end do
1147 : !$acc end parallel loop
1148 : end do
1149 :
1150 : ! A zero-flux boundary condition at the top of the model, d(xm)/dz = 0,
1151 : ! means that x'w' at the top model level is 0,
1152 : ! since x'w' = - K_zm * d(xm)/dz.
1153 : !$acc parallel loop gang vector collapse(2) default(present)
1154 8470656 : do j = 1, edsclr_dim
1155 135899856 : do i = 1, ngrdcol
1156 135546912 : wpedsclrp(i,nz,j) = zero
1157 : end do
1158 : end do
1159 : !$acc end parallel loop
1160 :
1161 : ! Compute the implicit portion of the xm (eddy-scalar) equations.
1162 : ! Build the left-hand side matrix.
1163 : call windm_edsclrm_lhs( nz, ngrdcol, gr, dt, & ! intent(in)
1164 : lhs_ma_zt, lhs_diff, & ! intent(in)
1165 : wind_speed, u_star_sqd, & ! intent(in)
1166 : rho_ds_zm, invrs_rho_ds_zt, & ! intent(in)
1167 : l_implemented, l_imp_sfc_momentum_flux, & ! intent(in)
1168 352944 : lhs ) ! intent(out)
1169 :
1170 : ! Decompose and back substitute for all eddy-scalar variables
1171 : call windm_edsclrm_solve( nz, ngrdcol, gr, edsclr_dim, 0, & ! intent(in)
1172 : tridiag_solve_method, & ! intent(in)
1173 : stats_metadata, & ! intent(in)
1174 : stats_sfc, & ! intent(inout)
1175 : lhs, rhs, & ! intent(inout)
1176 352944 : solution ) ! intent(out)
1177 :
1178 352944 : if ( clubb_at_least_debug_level( 0 ) ) then
1179 352944 : if ( err_code == clubb_fatal_error ) then
1180 0 : write(fstderr,*) "Fatal error solving for eddsclrm"
1181 : end if
1182 : end if
1183 :
1184 : !----------------------------------------------------------------
1185 : ! Update Eddy-diff. Passive Scalars
1186 : !----------------------------------------------------------------
1187 : !$acc parallel loop gang vector collapse(3) default(present)
1188 8470656 : do j = 1, edsclr_dim
1189 698476176 : do k = 1, nz
1190 11529605232 : do i = 1, ngrdcol
1191 11521487520 : edsclrm(i,k,j) = solution(i,k,j)
1192 : end do
1193 : end do
1194 : end do
1195 : !$acc end parallel loop
1196 :
1197 : ! The value of edsclrm(1) is located below the model surface and does not
1198 : ! effect the rest of the model. The value of edsclrm(1) is simply set to
1199 : ! the value of edsclrm(2) after the equation matrix has been solved.
1200 : !$acc parallel loop gang vector collapse(2) default(present)
1201 8470656 : do edsclr = 1, edsclr_dim
1202 135899856 : do i = 1, ngrdcol
1203 135546912 : edsclrm(i,1,edsclr) = edsclrm(i,2,edsclr)
1204 : end do
1205 : end do
1206 : !$acc end parallel loop
1207 :
1208 352944 : if ( l_lmm_stepping ) then
1209 : !$acc parallel loop gang vector collapse(3) default(present)
1210 0 : do j = 1, edsclr_dim
1211 0 : do k = 1, nz
1212 0 : do i = 1, ngrdcol
1213 0 : edsclrm(i,k,j) = one_half * ( edsclrm_old(i,k,j) + edsclrm(i,k,j) )
1214 : end do
1215 : end do
1216 : end do
1217 : !$acc end parallel loop
1218 : endif ! l_lmm_stepping
1219 :
1220 : ! Second part of momentum (implicit component)
1221 :
1222 : ! Solve for x'w' at all intermediate model levels.
1223 : ! A Crank-Nicholson timestep is used.
1224 8470656 : do edsclr = 1, edsclr_dim
1225 :
1226 : call calc_xpwp( nz, ngrdcol, gr, &
1227 : Kmh_zm, edsclrm(:,:,edsclr), &
1228 8117712 : xpwp )
1229 :
1230 : !$acc parallel loop gang vector collapse(2) default(present)
1231 682240752 : do k = 2, nz-1
1232 11258511408 : do i = 1, ngrdcol
1233 11250393696 : wpedsclrp(i,k,edsclr) = -one_half * xpwp(i,k)
1234 : end do
1235 : end do
1236 : !$acc end parallel loop
1237 : end do
1238 :
1239 : ! Note that the w'edsclr' terms are not clipped, since we don't compute
1240 : ! the variance of edsclr anywhere. -dschanen 7 Oct 2008
1241 :
1242 : endif
1243 :
1244 352944 : if ( clubb_at_least_debug_level( 0 ) ) then
1245 352944 : if ( err_code == clubb_fatal_error ) then
1246 :
1247 : !$acc update host( wm_zt, Km_zm, ug, vg, um_ref, vm_ref, wp2, &
1248 : !$acc up2, vp2, um_forcing, vm_forcing, edsclrm_forcing, &
1249 : !$acc fcor, um_old, um, vm_old, vm, edsclrm_old, &
1250 : !$acc edsclrm, upwp, vpwp, wpedsclrp )
1251 :
1252 0 : write(fstderr,*) "Error in advance_windm_edsclrm"
1253 :
1254 0 : write(fstderr,*) "intent(in)"
1255 :
1256 0 : write(fstderr,*) "dt = ", dt
1257 0 : write(fstderr,*) "wm_zt = ", wm_zt
1258 0 : write(fstderr,*) "Km_zm = ", Km_zm
1259 0 : write(fstderr,*) "ug = ", ug
1260 0 : write(fstderr,*) "vg = ", vg
1261 0 : write(fstderr,*) "um_ref = ", um_ref
1262 0 : write(fstderr,*) "vm_ref = ", vm_ref
1263 0 : write(fstderr,*) "wp2 = ", wp2
1264 0 : write(fstderr,*) "up2 = ", up2
1265 0 : write(fstderr,*) "vp2 = ", vp2
1266 0 : write(fstderr,*) "um_forcing = ", um_forcing
1267 0 : write(fstderr,*) "vm_forcing = ", vm_forcing
1268 0 : do edsclr = 1, edsclr_dim
1269 0 : write(fstderr,*) "edsclrm_forcing # = ",edsclr, edsclrm_forcing
1270 : end do
1271 0 : write(fstderr,*) "fcor = ", fcor
1272 0 : write(fstderr,*) "l_implemented = ", l_implemented
1273 :
1274 0 : write(fstderr,*) "intent(inout)"
1275 :
1276 0 : if ( l_lmm_stepping ) &
1277 0 : write(fstderr,*) "um (pre-solve) = ", um_old
1278 0 : write(fstderr,*) "um = ", um
1279 0 : if ( l_lmm_stepping ) &
1280 0 : write(fstderr,*) "vm (pre-solve) = ", vm_old
1281 0 : write(fstderr,*) "vm = ", vm
1282 0 : do edsclr = 1, edsclr_dim
1283 0 : if ( l_lmm_stepping ) &
1284 0 : write(fstderr,*) "edsclrm (pre-solve) # ", edsclr, "=", edsclrm_old(:,:,edsclr)
1285 0 : write(fstderr,*) "edsclrm # ", edsclr, "=", edsclrm(:,:,edsclr)
1286 : end do
1287 0 : write(fstderr,*) "upwp = ", upwp
1288 0 : write(fstderr,*) "vpwp = ", vpwp
1289 0 : write(fstderr,*) "wpedsclrp = ", wpedsclrp
1290 :
1291 0 : return
1292 : end if
1293 : end if
1294 :
1295 : !$acc exit data delete( um_old, vm_old, um_tndcy, vm_tndcy, &
1296 : !$acc upwp_chnge, vpwp_chnge, lhs, rhs, solution, wind_speed, &
1297 : !$acc wind_speed_pert, u_star_sqd, u_star_sqd_pert, &
1298 : !$acc nu_zero, lhs_diff, lhs_ma_zt, Km_zt, Kmh_zt, &
1299 : !$acc Km_zm_p_nu10, xpwp )
1300 :
1301 : !$acc exit data if( edsclr_dim > 0) delete( edsclrm_old )
1302 :
1303 : return
1304 :
1305 : end subroutine advance_windm_edsclrm
1306 :
1307 : !=============================================================================
1308 352944 : subroutine windm_edsclrm_solve( nz, ngrdcol, gr, nrhs, ixm_matrix_condt_num, &
1309 : tridiag_solve_method, &
1310 : stats_metadata, &
1311 352944 : stats_sfc, &
1312 352944 : lhs, rhs, solution )
1313 :
1314 : ! Note: In the "Description" section of this subroutine, the variable
1315 : ! "invrs_dzm" will be written as simply "dzm", and the variable
1316 : ! "invrs_dzt" will be written as simply "dzt". This is being done as
1317 : ! as device to save space and to make some parts of the description
1318 : ! more readable. This change does not pertain to the actual code.
1319 :
1320 : ! Description:
1321 : ! Solves the horizontal wind or eddy-scalar time-tendency equation, and
1322 : ! diagnoses the turbulent flux. A Crank-Nicholson time-stepping algorithm
1323 : ! is used in solving the turbulent advection term and in diagnosing the
1324 : ! turbulent flux.
1325 : !
1326 : ! The rate of change of an eddy-scalar quantity, xm, is:
1327 : !
1328 : ! d(xm)/dt = - w * d(xm)/dz - (1/rho_ds) * d( rho_ds * x'w' )/dz
1329 : ! + xm_forcings.
1330 : !
1331 : !
1332 : ! The Turbulent Advection Term
1333 : ! ----------------------------
1334 : !
1335 : ! The above equation contains a turbulent advection term:
1336 : !
1337 : ! - (1/rho_ds) * d( rho_ds * x'w' )/dz;
1338 : !
1339 : ! where the momentum flux, x'w', is closed using a down gradient approach:
1340 : !
1341 : ! x'w' = - K_zm * d(xm)/dz.
1342 : !
1343 : ! The turbulent advection term becomes:
1344 : !
1345 : ! + (1/rho_ds) * d [ rho_ds * K_zm * d(xm)/dz ] / dz;
1346 : !
1347 : ! which is the same as a standard eddy-diffusion term (if "rho_ds * K_zm" in
1348 : ! the term above is substituted for "K_zm" in a standard eddy-diffusion
1349 : ! term, and if the standard eddy-diffusion term is multiplied by
1350 : ! "1/rho_ds"). Thus, the turbulent advection term is treated and solved in
1351 : ! the same way that a standard eddy-diffusion term would be solved. The
1352 : ! term is discretized as follows:
1353 : !
1354 : ! The values of xm are found on the thermodynamic levels, while the values
1355 : ! of K_zm are found on the momentum levels. Additionally, the values of
1356 : ! rho_ds_zm are found on the momentum levels, and the values of
1357 : ! invrs_rho_ds_zt are found on the thermodynamic levels. The
1358 : ! derivatives (d/dz) of xm are taken over the intermediate momentum levels.
1359 : ! At the intermediate momentum levels, d(xm)/dz is multiplied by K_zm and by
1360 : ! rho_ds_zm. Then, the derivative of the whole mathematical expression is
1361 : ! taken over the central thermodynamic level, where it is multiplied by
1362 : ! invrs_rho_ds_zt, which yields the desired result.
1363 : !
1364 : ! ---xm(kp1)----------------------------------------------------- t(k+1)
1365 : !
1366 : ! ===========d(xm)/dz===K_zm(k)=====rho_ds_zm(k)================= m(k)
1367 : !
1368 : ! ---xm(k)---invrs_rho_ds_zt---d[rho_ds_zm*K_zm*d(xm)/dz]/dz----- t(k)
1369 : !
1370 : ! ===========d(xm)/dz===K_zm(km1)===rho_ds_zm(km1)=============== m(k-1)
1371 : !
1372 : ! ---xm(km1)----------------------------------------------------- t(k-1)
1373 : !
1374 : ! The vertical indices t(k+1), m(k), t(k), m(k-1), and t(k-1) correspond
1375 : ! with altitudes zt(k+1), zm(k), zt(k), zm(k-1), and zt(k-1), respectively.
1376 : ! The letter "t" is used for thermodynamic levels and the letter "m" is used
1377 : ! for momentum levels.
1378 : !
1379 : ! dzt(k) = 1 / ( zm(k) - zm(k-1) )
1380 : ! dzm(k) = 1 / ( zt(k+1) - zt(k) )
1381 : ! dzm(k-1) = 1 / ( zt(k) - zt(k-1) )
1382 : !
1383 : ! The vertically discretized form of the turbulent advection term (treated
1384 : ! as an eddy diffusion term) is written out as:
1385 : !
1386 : ! + invrs_rho_ds_zt(k)
1387 : ! * dzt(k)
1388 : ! * [ rho_ds_zm(k) * K_zm(k) * dzm(k) * ( xm(k+1) - xm(k) )
1389 : ! - rho_ds_zm(k-1) * K_zm(k-1) * dzm(k-1) * ( xm(k) - xm(k-1) ) ].
1390 : !
1391 : ! For this equation, a Crank-Nicholson (semi-implicit) diffusion scheme is
1392 : ! used to solve the (1/rho_ds) * d [ rho_ds * K_zm * d(xm)/dz ] / dz
1393 : ! eddy-diffusion term. The discretized implicit portion of the term is
1394 : ! written out as:
1395 : !
1396 : ! + (1/2) * invrs_rho_ds_zt(k)
1397 : ! * dzt(k)
1398 : ! * [ rho_ds_zm(k) * K_zm(k)
1399 : ! * dzm(k) * ( xm(k+1,<t+1>) - xm(k,<t+1>) )
1400 : ! - rho_ds_zm(k-1) * K_zm(k-1)
1401 : ! * dzm(k-1) * ( xm(k,<t+1>) - xm(k-1,<t+1>) ) ].
1402 : !
1403 : ! Note: When the implicit term is brought over to the left-hand side,
1404 : ! the sign is reversed and the leading "+" in front of the term
1405 : ! is changed to a "-".
1406 : !
1407 : ! The discretized explicit portion of the term is written out as:
1408 : !
1409 : ! + (1/2) * invrs_rho_ds_zt(k)
1410 : ! * dzt(k)
1411 : ! * [ rho_ds_zm(k) * K_zm(k)
1412 : ! * dzm(k) * ( xm(k+1,<t>) - xm(k,<t>) )
1413 : ! - rho_ds_zm(k-1) * K_zm(k-1)
1414 : ! * dzm(k-1) * ( xm(k,<t>) - xm(k-1,<t>) ) ].
1415 : !
1416 : ! Timestep index (t) stands for the index of the current timestep, while
1417 : ! timestep index (t+1) stands for the index of the next timestep, which is
1418 : ! being advanced to in solving the d(xm)/dt equation.
1419 : !
1420 : !
1421 : ! Boundary Conditions:
1422 : !
1423 : ! An eddy-scalar quantity is not allowed to flux out the upper boundary.
1424 : ! Thus, a zero-flux boundary condition is used for the upper boundary in the
1425 : ! eddy-diffusion equation.
1426 : !
1427 : ! The lower boundary condition is much more complicated. It is neither a
1428 : ! zero-flux nor a fixed-point boundary condition. Rather, it is a
1429 : ! fixed-flux boundary condition. This term is a turbulent advection term,
1430 : ! but with the eddy-scalars, the only value of x'w' relevant in solving the
1431 : ! d(xm)/dt equation is the value of x'w' at the surface (the first momentum
1432 : ! level), which is written as x'w'|_sfc.
1433 : !
1434 : ! 1) x'w' surface flux; generalized explicit form
1435 : !
1436 : ! The x'w' surface flux is applied to the d(xm)/dt equation through the
1437 : ! turbulent advection term, which is:
1438 : !
1439 : ! - (1/rho_ds) * d( rho_ds * x'w' )/dz.
1440 : !
1441 : ! At most vertical levels, a substitution can be made for x'w', such
1442 : ! that:
1443 : !
1444 : ! x'w' = - K_zm * d(xm)/dz.
1445 : !
1446 : ! However, the same substitution cannot be made at the surface (momentum
1447 : ! level 1), as x'w'|_sfc is a surface flux that is explicitly computed
1448 : ! elsewhere in the model code.
1449 : !
1450 : ! The lower boundary condition, which in this case needs to be applied to
1451 : ! the d(xm)/dt equation at level 2, is discretized as follows:
1452 : !
1453 : ! --xm(3)------------------------------------------------------- t(3)
1454 : !
1455 : ! ========[x'w'(2) = -K_zm(2)*d(xm)/dz]===rho_ds_zm(2)========== m(2)
1456 : !
1457 : ! --xm(2)---invrs_rho_ds_zt(2)---d[rho_ds_zm*K_zm*d(xm)/dz]/dz-- t(2)
1458 : !
1459 : ! ========[x'w'|_sfc]=====================rho_ds_zm(1)========== m(1) sfc
1460 : !
1461 : ! --xm(1)-------(below surface; not applicable)----------------- t(1)
1462 : !
1463 : ! where "sfc" is the level of the model surface or lower boundary.
1464 : !
1465 : ! The vertically discretized form of the turbulent advection term
1466 : ! (treated as an eddy diffusion term), with the explicit surface flux,
1467 : ! x'w'|_sfc, in place, is written out as:
1468 : !
1469 : ! - invrs_rho_ds_zt(2)
1470 : ! * dzt(2) * [ rho_ds_zm(2) * x'w'(2) - rho_ds_zm(1) * x'w'|_sfc ];
1471 : !
1472 : ! which can be re-written as:
1473 : !
1474 : ! + invrs_rho_ds_zt(2)
1475 : ! * dzt(2)
1476 : ! * [ rho_ds_zm(2) * K_zm(2) * dzm(2) * ( xm(3) - xm(2) )
1477 : ! + rho_ds_zm(1) * x'w'|_sfc ];
1478 : !
1479 : ! which can be re-written again as:
1480 : !
1481 : ! + invrs_rho_ds_zt(2)
1482 : ! * dzt(2)
1483 : ! * rho_ds_zm(2) * K_zm(2) * dzm(2) * ( xm(3) - xm(2) )
1484 : ! + invrs_rho_ds_zt(2)
1485 : ! * dzt(2)
1486 : ! * rho_ds_zm(1) * x'w'|_sfc.
1487 : !
1488 : ! For this equation, a Crank-Nicholson (semi-implicit) diffusion scheme
1489 : ! is used to solve the (1/rho_ds) * d [ rho_ds * K_zm * d(xm)/dz ] / dz
1490 : ! eddy-diffusion term. The discretized implicit portion of the term is
1491 : ! written out as:
1492 : !
1493 : ! + (1/2) * invrs_rho_ds_zt(2)
1494 : ! * dzt(2)
1495 : ! * [ rho_ds_zm(2) * K_zm(2)
1496 : ! * dzm(2) * ( xm(3,<t+1>) - xm(2,<t+1>) ) ].
1497 : !
1498 : ! Note: When the implicit term is brought over to the left-hand side,
1499 : ! the sign is reversed and the leading "+" in front of the term
1500 : ! is changed to a "-".
1501 : !
1502 : ! The discretized explicit portion of the term is written out as:
1503 : !
1504 : ! + (1/2) * invrs_rho_ds_zt(2)
1505 : ! * dzt(2)
1506 : ! * [ rho_ds_zm(2) * K_zm(2)
1507 : ! * dzm(2) * ( xm(3,<t>) - xm(2,<t>) ) ]
1508 : ! + invrs_rho_ds_zt(2)
1509 : ! * dzt(2)
1510 : ! * rho_ds_zm(1) * x'w'|_sfc.
1511 : !
1512 : ! Note: The x'w'|_sfc portion of the term written above has been pulled
1513 : ! away from the rest of the explicit form written above because
1514 : ! the (1/2) factor due to Crank-Nicholson time_stepping does not
1515 : ! apply to it, as there isn't an implicit portion for x'w'|_sfc.
1516 : !
1517 : ! Timestep index (t) stands for the index of the current timestep, while
1518 : ! timestep index (t+1) stands for the index of the next timestep, which
1519 : ! is being advanced to in solving the d(xm)/dt equation.
1520 : !
1521 : ! 2) x'w' surface flux; implicit form for momentum fluxes u'w' and v'w'
1522 : !
1523 : ! The x'w' surface flux is applied to the d(xm)/dt equation through the
1524 : ! turbulent advection term, which is:
1525 : !
1526 : ! - (1/rho_ds) * d( rho_ds * x'w' )/dz.
1527 : !
1528 : ! At most vertical levels, a substitution can be made for x'w', such
1529 : ! that:
1530 : !
1531 : ! x'w' = - K_zm * d(xm)/dz.
1532 : !
1533 : ! However, the same substitution cannot be made at the surface (momentum
1534 : ! level 1), as x'w'|_sfc is a surface momentum flux that is found by the
1535 : ! following equation:
1536 : !
1537 : ! x'w'|_sfc = - [ u_star^2 / sqrt( um^2 + vm^2 ) ] * xm;
1538 : !
1539 : ! where x'w'|_sfc and xm are either u'w'|_sfc and um, respectively, or
1540 : ! v'w'|_sfc and vm, respectively (um and vm are located at the first
1541 : ! thermodynamic level above the surface, which is thermodynamic level 2),
1542 : ! sqrt( um^2 + vm^2 ) is the wind speed (also at thermodynamic level 2),
1543 : ! and u_star is defined as:
1544 : !
1545 : ! u_star = ( u'w'|_sfc^2 + v'w'|_sfc^2 )^(1/4);
1546 : !
1547 : ! and thus u_star^2 is defined as:
1548 : !
1549 : ! u_star^2 = sqrt( u'w'|_sfc^2 + v'w'|_sfc^2 ).
1550 : !
1551 : ! The value of u_star is either set to a constant value or computed
1552 : ! (through function diag_ustar) based on the surface wind speed, the
1553 : ! height above surface of the surface wind speed (as compared to the
1554 : ! roughness height), and the buoyancy flux at the surface. Either way,
1555 : ! u_star is computed elsewhere in the model, and the values of u'w'|_sfc
1556 : ! and v'w'|_sfc are based on it and computed along with it. The values
1557 : ! of u'w'|_sfc and v'w'|_sfc are then passed into advance_clubb_core,
1558 : ! and are eventually passed into advance_windm_edsclrm. In subroutine
1559 : ! advance_windm_edsclrm, the value of u_star_sqd is then recomputed
1560 : ! based on u'w'|_sfc and v'w'|_sfc. The value of sqrt( u_star_sqd ) is
1561 : ! consistent with the value of the original computation of u_star.
1562 : !
1563 : ! The equation listed above is substituted for x'w'|_sfc. The lower
1564 : ! boundary condition, which in this case needs to be applied to the
1565 : ! d(xm)/dt equation at level 2, is discretized as follows:
1566 : !
1567 : ! --xm(3)------------------------------------------------------- t(3)
1568 : !
1569 : ! ===[x'w'(2) = -K_zm(2)*d(xm)/dz]=================rho_ds_zm(2)= m(2)
1570 : !
1571 : ! --xm(2)---invrs_rho_ds_zt(2)---d[rho_ds_zm*K_zm*d(xm)/dz]/dz-- t(2)
1572 : !
1573 : ! ===[x'w'|_sfc = -[u_star^2/sqrt(um^2+vm^2)]*xm]==rho_ds_zm(1)= m(1) sfc
1574 : !
1575 : ! --xm(1)-------(below surface; not applicable)----------------- t(1)
1576 : !
1577 : ! where "sfc" is the level of the model surface or lower boundary.
1578 : !
1579 : ! The vertically discretized form of the turbulent advection term
1580 : ! (treated as an eddy diffusion term), with the implicit surface momentum
1581 : ! flux in place, is written out as:
1582 : !
1583 : ! - invrs_rho_ds_zt(2)
1584 : ! * dzt(2) * [ rho_ds_zm(2) * x'w'(2) - rho_ds_zm(1) * x'w'|_sfc ];
1585 : !
1586 : ! which can be re-written as:
1587 : !
1588 : ! - invrs_rho_ds_zt(2)
1589 : ! * dzt(2)
1590 : ! * [ rho_ds_zm(2)
1591 : ! * { - K_zm(2) * dzm(2) * ( xm(3) - xm(2) ) }
1592 : ! - rho_ds_zm(1)
1593 : ! * { - [ u_star^2 / sqrt( um(2)^2 + vm(2)^2 ) ] * xm(2) } ];
1594 : !
1595 : ! which can be re-written as:
1596 : !
1597 : ! + invrs_rho_ds_zt(2)
1598 : ! * dzt(2)
1599 : ! * rho_ds_zm(2) * K_zm(2) * dzm(2) * ( xm(3) - xm(2) )
1600 : ! - invrs_rho_ds_zt(2)
1601 : ! * dzt(2)
1602 : ! * rho_ds_zm(1) * [ u_star^2 / sqrt( um(2)^2 + vm(2)^2 ) ] * xm(2).
1603 : !
1604 : ! For this equation, a Crank-Nicholson (semi-implicit) diffusion scheme
1605 : ! is used to solve the (1/rho_ds) * d [ rho_ds * K_zm * d(xm)/dz ] / dz
1606 : ! eddy-diffusion term. The discretized implicit portion of the term is
1607 : ! written out as:
1608 : !
1609 : ! + (1/2) * invrs_rho_ds_zt(2)
1610 : ! * dzt(2)
1611 : ! * [ rho_ds_zm(2) * K_zm(2)
1612 : ! * dzm(2) * ( xm(3,<t+1>) - xm(2,<t+1>) ) ]
1613 : ! - invrs_rho_ds_zt(2)
1614 : ! * dzt(2)
1615 : ! * rho_ds_zm(1)
1616 : ! * [u_star^2/sqrt( um(2,<t>)^2 + vm(2,<t>)^2 )] * xm(2,<t+1>).
1617 : !
1618 : ! Note: When the implicit term is brought over to the left-hand side,
1619 : ! the signs are reversed and the leading "+" in front of the first
1620 : ! part of the term is changed to a "-", while the leading "-" in
1621 : ! front of the second part of the term is changed to a "+".
1622 : !
1623 : ! Note: The x'w'|_sfc portion of the term written above has been pulled
1624 : ! away from the rest of the implicit form written above because
1625 : ! the (1/2) factor due to Crank-Nicholson time_stepping does not
1626 : ! apply to it. The x'w'|_sfc portion of the term is treated
1627 : ! completely implicitly in order to enhance numerical stability.
1628 : !
1629 : ! The discretized explicit portion of the term is written out as:
1630 : !
1631 : ! + (1/2) * invrs_rho_ds_zt(2)
1632 : ! * dzt(2)
1633 : ! * [ rho_ds_zm(2) * K_zm(2)
1634 : ! * dzm(2) * ( xm(3,<t>) - xm(2,<t>) ) ].
1635 : !
1636 : ! Timestep index (t) stands for the index of the current timestep, while
1637 : ! timestep index (t+1) stands for the index of the next timestep, which
1638 : ! is being advanced to in solving the d(xm)/dt equation.
1639 : !
1640 : !
1641 : ! The lower boundary condition for the implicit and explicit portions of the
1642 : ! turbulent advection term, without the x'w'|_sfc portion of the term, can
1643 : ! easily be invoked by using the zero-flux boundary conditions found in the
1644 : ! generalized diffusion function (function diffusion_zt_lhs), which is used
1645 : ! for many other equations in this model. Either the generalized explicit
1646 : ! surface flux needs to be added onto the explicit term after the diffusion
1647 : ! function has been called from subroutine windm_edsclrm_rhs, or the
1648 : ! implicit momentum surface flux needs to be added onto the implicit term
1649 : ! after the diffusion function has been called from subroutine
1650 : ! windm_edsclrm_lhs. However, all other equations in this model that use
1651 : ! zero-flux diffusion have level 1 as the level to which the lower boundary
1652 : ! condition needs to be applied. Thus, an adjuster will have to be used at
1653 : ! level 2 to call diffusion_zt_lhs with level 1 as the input level (the last
1654 : ! variable being passed in during the function call). However, the other
1655 : ! variables passed in (rho_ds_zm*K_zm, gr%dzt(1,:), and gr%dzm(1,:) variables) will
1656 : ! have to be passed in as solving for level 2.
1657 : !
1658 : ! The value of xm(1) is located below the model surface and does not effect
1659 : ! the rest of the model. Since xm can be either a horizontal wind component
1660 : ! or a generic eddy scalar quantity, the value of xm(1) is simply set to the
1661 : ! value of xm(2) after the equation matrix has been solved.
1662 : !
1663 : !
1664 : ! Conservation Properties:
1665 : !
1666 : ! When a fixed-flux lower boundary condition is used (combined with a
1667 : ! zero-flux upper boundary condition), this technique of discretizing the
1668 : ! turbulent advection term (treated as an eddy-diffusion term) leads to
1669 : ! conservative differencing. When the implicit momentum surface flux is
1670 : ! either zero or not used, the column totals for each column in the
1671 : ! left-hand side matrix (for the turbulent advection term) should be equal
1672 : ! to 0. Otherwise, the column total for the second column will be equal to
1673 : ! rho_ds_zm(1) * x'w'|_sfc<t+1>. When the generalized explicit surface
1674 : ! flux is either zero or not used, the column total for the right-hand side
1675 : ! vector (for the turbulent advection term) should be equal to 0.
1676 : ! Otherwise, the column total for the right-hand side vector (for the
1677 : ! turbulent advection term) will be equal to rho_ds_zm(1) * x'w'|_sfc<t>.
1678 : ! This ensures that the total amount of quantity xm over the entire vertical
1679 : ! domain is only changed by the surface flux (neglecting any forcing terms).
1680 : ! The total amount of change is equal to rho_ds_zm(1) * x'w'|_sfc.
1681 : !
1682 : ! To see that this conservation law is satisfied by the left-hand side
1683 : ! matrix, compute the turbulent advection (treated as eddy diffusion) of xm,
1684 : ! neglecting any implicit momentum surface flux, multiply by rho_ds_zt, and
1685 : ! integrate vertically. In discretized matrix notation (where "i" stands
1686 : ! for the matrix column and "j" stands for the matrix row):
1687 : !
1688 : ! 0 = Sum_j Sum_i
1689 : ! (rho_ds_zt)_i ( 1/dzt )_i
1690 : ! ( 0.5_core_rknd * (1/rho_ds_zt) * dzt * (rho_ds_zm*K_zm*dzm) )_ij (xm<t+1>)_j.
1691 : !
1692 : ! The left-hand side matrix,
1693 : ! ( 0.5_core_rknd * (1/rho_ds_zt) * dzt * (rho_ds_zm*K_zm*dzm) )_ij, is partially
1694 : ! written below. The sum over i in the above equation removes (1/rho_ds_zt)
1695 : ! and dzt everywhere from the matrix below. The sum over j leaves the
1696 : ! column totals that are desired, which are 0.
1697 : !
1698 : ! Left-hand side matrix contributions from the turbulent advection term
1699 : ! (treated as an eddy-diffusion term using a Crank-Nicholson timestep);
1700 : ! first five vertical levels:
1701 : !
1702 : ! ------------------------------------------------------------------------------->
1703 : !k=1 | 0 0 0 0
1704 : ! |
1705 : !k=2 | 0 +0.5* -0.5* 0
1706 : ! | (1/rho_ds_zt(k))* (1/rho_ds_zt(k))*
1707 : ! | dzt(k)* dzt(k)*
1708 : ! | rho_ds_zm(k)* rho_ds_zm(k)*
1709 : ! | K_zm(k)*dzm(k) K_zm(k)*dzm(k)
1710 : ! |
1711 : !k=3 | 0 -0.5* +0.5* -0.5*
1712 : ! | (1/rho_ds_zt(k))* (1/rho_ds_zt(k))* (1/rho_ds_zt(k))*
1713 : ! | dzt(k)* dzt(k)* dzt(k)*
1714 : ! | rho_ds_zm(k-1)* [ rho_ds_zm(k)* rho_ds_zm(k)*
1715 : ! | K_zm(k-1)*dzm(k-1) K_zm(k)*dzm(k) K_zm(k)*dzm(k)
1716 : ! | +rho_ds_zm(k-1)*
1717 : ! | K_zm(k-1)*dzm(k-1) ]
1718 : ! |
1719 : !k=4 | 0 0 -0.5* +0.5*
1720 : ! | (1/rho_ds_zt(k))* (1/rho_ds_zt(k))*
1721 : ! | dzt(k)* dzt(k)*
1722 : ! | rho_ds_zm(k-1)* [ rho_ds_zm(k)*
1723 : ! | K_zm(k-1)*dzm(k-1) K_zm(k)*dzm(k)
1724 : ! | +rho_ds_zm(k-1)*
1725 : ! | K_zm(k-1)*dzm(k-1) ]
1726 : ! |
1727 : !k=5 | 0 0 0 -0.5*
1728 : ! | (1/rho_ds_zt(k))*
1729 : ! | dzt(k)*
1730 : ! | rho_ds_zm(k-1)*
1731 : ! | K_zm(k-1)*dzm(k-1)
1732 : ! \ /
1733 : !
1734 : ! Note: The superdiagonal term from level 4 and both the main diagonal and
1735 : ! superdiagonal terms from level 5 are not shown on this diagram.
1736 : !
1737 : ! Note: If an implicit momentum surface flux is used, an additional term,
1738 : ! + (1/rho_ds_zt(2)) * dzt(2) * rho_ds_zm(1)
1739 : ! * [ u_star^2 / sqrt( um(2,<t>)^2 + vm(2,<t>)^2 ) ], is added to
1740 : ! row 2 (k=2), column 2.
1741 : !
1742 : ! To see that the above conservation law is satisfied by the right-hand side
1743 : ! vector, compute the turbulent advection (treated as eddy diffusion) of xm,
1744 : ! neglecting any generalized explicit surface flux, multiply by rho_ds_zt,
1745 : ! and integrate vertically. In discretized matrix notation (where "i"
1746 : ! stands for the matrix column and "j" stands for the matrix row):
1747 : !
1748 : ! 0 = Sum_j Sum_i (rho_ds_zt)_i ( 1/dzt )_i ( rhs_vector )_j.
1749 : !
1750 : ! The right-hand side vector, ( rhs_vector )_j, is partially written below.
1751 : ! The sum over i in the above equation removes (1/rho_ds_zt) and dzt
1752 : ! everywhere from the vector below. The sum over j leaves the column total
1753 : ! that is desired, which is 0.
1754 : !
1755 : ! Right-hand side vector contributions from the turbulent advection term
1756 : ! (treated as an eddy-diffusion term using a Crank-Nicholson timestep);
1757 : ! first five vertical levels:
1758 : !
1759 : ! --------------------------------------------
1760 : !k=1 | 0 |
1761 : ! | |
1762 : ! | |
1763 : !k=2 | +0.5*(1/rho_ds_zt(k))* |
1764 : ! | dzt(k)* |
1765 : ! | [ rho_ds_zm(k)*K_zm(k)* |
1766 : ! | dzm(k)*(xm(k+1,<t>)-xm(k,<t>)) ] |
1767 : ! | |
1768 : !k=3 | +0.5*(1/rho_ds_zt(k))* |
1769 : ! | dzt(k)* |
1770 : ! | [ rho_ds_zm(k)*K_zm(k)* |
1771 : ! | dzm(k)*(xm(k+1,<t>)-xm(k,<t>)) |
1772 : ! | -rho_ds_zm(k-1)*K_zm(k-1)* |
1773 : ! | dzm(k-1)*(xm(k,<t>)-xm(k-1,<t>)) ] |
1774 : ! | |
1775 : !k=4 | +0.5*(1/rho_ds_zt(k))* |
1776 : ! | dzt(k)* |
1777 : ! | [ rho_ds_zm(k)*K_zm(k)* |
1778 : ! | dzm(k)*(xm(k+1,<t>)-xm(k,<t>)) |
1779 : ! | -rho_ds_zm(k-1)*K_zm(k-1)* |
1780 : ! | dzm(k-1)*(xm(k,<t>)-xm(k-1,<t>)) ] |
1781 : ! | |
1782 : !k=5 | +0.5*(1/rho_ds_zt(k))* |
1783 : ! | dzt(k)* |
1784 : ! | [ rho_ds_zm(k)*K_zm(k)* |
1785 : ! | dzm(k)*(xm(k+1,<t>)-xm(k,<t>)) |
1786 : ! | -rho_ds_zm(k-1)*K_zm(k-1)* |
1787 : ! | dzm(k-1)*(xm(k,<t>)-xm(k-1,<t>)) ] |
1788 : ! \ / \ /
1789 : !
1790 : ! Note: If a generalized explicit surface flux is used, an additional term,
1791 : ! + (1/rho_ds_zt(2)) * dzt(2) * rho_ds_zm(1) * x'w'|_sfc, is added to
1792 : ! row 2 (k=2).
1793 : !
1794 : ! Note: Only the contributions by the turbulent advection term are shown
1795 : ! for both the left-hand side matrix and the right-hand side vector.
1796 : ! There are more terms in the equation, and thus more factors to be
1797 : ! added to both the left-hand side matrix (such as time tendency and
1798 : ! mean advection) and the right-hand side vector (such as xm
1799 : ! forcings). The left-hand side matrix is set-up so that a singular
1800 : ! matrix is not encountered.
1801 :
1802 : ! References:
1803 : ! Eqn. 8 & 9 on p. 3545 of
1804 : ! ``A PDF-Based Model for Boundary Layer Clouds. Part I:
1805 : ! Method and Model Description'' Golaz, et al. (2002)
1806 : ! JAS, Vol. 59, pp. 3540--3551.
1807 : !-----------------------------------------------------------------------
1808 :
1809 : use grid_class, only: &
1810 : grid ! Type
1811 :
1812 : use matrix_solver_wrapper, only: &
1813 : tridiag_solve ! Procedure(s)
1814 :
1815 : use stats_variables, only: &
1816 : stats_metadata_type
1817 :
1818 : use stats_type_utilities, only: &
1819 : stat_update_var_pt ! Subroutine
1820 :
1821 : use clubb_precision, only: &
1822 : core_rknd ! Variable(s)
1823 :
1824 : use stats_type, only: stats ! Type
1825 :
1826 : implicit none
1827 :
1828 : ! Constant parameters
1829 :
1830 : integer, parameter :: &
1831 : kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
1832 : k_tdiag = 2, & ! Thermodynamic main diagonal index.
1833 : km1_tdiag = 3 ! Thermodynamic subdiagonal index.
1834 :
1835 : ! ------------------------ Input Variables ------------------------
1836 : integer, intent(in) :: &
1837 : nz, &
1838 : ngrdcol
1839 :
1840 : type (grid), target, intent(in) :: &
1841 : gr
1842 :
1843 : integer, intent(in) :: &
1844 : nrhs ! Number of right-hand side (explicit) vectors & Number of solution vectors.
1845 :
1846 : integer, intent(in) :: &
1847 : ixm_matrix_condt_num ! Stats index of the condition numbers
1848 :
1849 : integer, intent(in) :: &
1850 : tridiag_solve_method ! Specifier for method to solve tridiagonal systems
1851 :
1852 : type (stats_metadata_type), intent(in) :: &
1853 : stats_metadata
1854 :
1855 : ! ------------------------ Inout variables ------------------------
1856 : type (stats), target, dimension(ngrdcol), intent(inout) :: &
1857 : stats_sfc
1858 :
1859 : real( kind = core_rknd ), dimension(3,ngrdcol,nz), intent(inout) :: &
1860 : lhs ! Implicit contributions to um, vm, and eddy scalars [units vary]
1861 :
1862 : real( kind = core_rknd ), dimension(ngrdcol,nz,nrhs), intent(inout) :: &
1863 : rhs ! Right-hand side (explicit) contributions.
1864 :
1865 : ! ------------------------ Output variables ------------------------
1866 : real( kind = core_rknd ), dimension(ngrdcol,nz,nrhs), intent(out) :: &
1867 : solution ! Solution to the system of equations [units vary]
1868 :
1869 : ! ------------------------ Local variables ------------------------
1870 : real( kind = core_rknd ), dimension(ngrdcol) :: &
1871 705888 : rcond ! Estimate of the reciprocal of the condition number on the LHS matrix
1872 :
1873 : integer :: i
1874 :
1875 : ! ------------------------ Begin Code ------------------------
1876 :
1877 : ! Solve tridiagonal system for xm.
1878 352944 : if ( stats_metadata%l_stats_samp .and. ixm_matrix_condt_num > 0 ) then
1879 :
1880 : call tridiag_solve( "windm_edsclrm", tridiag_solve_method, & ! Intent(in)
1881 : ngrdcol, nz, nrhs, & ! Intent(in)
1882 : lhs, rhs, & ! Intent(inout)
1883 0 : solution, rcond ) ! Intent(out)
1884 :
1885 : ! Est. of the condition number of the variance LHS matrix
1886 0 : do i = 1, ngrdcol
1887 0 : call stat_update_var_pt( ixm_matrix_condt_num, 1, 1.0_core_rknd/rcond(i), & ! intent(in)
1888 0 : stats_sfc(i) ) ! intent(inout)
1889 : end do
1890 : else
1891 :
1892 : call tridiag_solve( "windm_edsclrm", tridiag_solve_method, & ! Intent(in)
1893 : ngrdcol, nz, nrhs, & ! Intent(in)
1894 : lhs, rhs, & ! Intent(inout)
1895 352944 : solution ) ! Intent(out)
1896 : end if
1897 :
1898 352944 : return
1899 : end subroutine windm_edsclrm_solve
1900 :
1901 : !=============================================================================
1902 0 : subroutine windm_edsclrm_implicit_stats( nz, solve_type, xm, &
1903 0 : invrs_dzt, &
1904 0 : lhs_diff, lhs_ma_zt, &
1905 0 : invrs_rho_ds_zt, u_star_sqd,&
1906 0 : rho_ds_zm, wind_speed, &
1907 : l_imp_sfc_momentum_flux, &
1908 : stats_metadata, &
1909 : stats_zt )
1910 :
1911 : ! Description:
1912 : ! Compute implicit contributions to um and vm
1913 :
1914 : ! References:
1915 : ! None
1916 : !-----------------------------------------------------------------------
1917 :
1918 : use constants_clubb, only: &
1919 : zero
1920 :
1921 : use stats_type_utilities, only: &
1922 : stat_end_update_pt, & ! Subroutines
1923 : stat_update_var_pt
1924 :
1925 : use clubb_precision, only: &
1926 : core_rknd
1927 :
1928 : use stats_type, only: &
1929 : stats ! Type
1930 :
1931 : use stats_variables, only: &
1932 : stats_metadata_type
1933 :
1934 : implicit none
1935 :
1936 : !---------------------- Input Variables ----------------------
1937 : integer, intent(in) :: &
1938 : nz
1939 :
1940 : integer, intent(in) :: &
1941 : solve_type ! Desc. of what is being solved for
1942 :
1943 : real( kind = core_rknd ), dimension(nz), intent(in) :: &
1944 : xm, & ! Computed value um or vm at <t+1> [m/s]
1945 : invrs_dzt ! The inverse spacing between momentum grid levels;
1946 : ! centered over thermodynamic grid levels.
1947 :
1948 : real( kind = core_rknd ), dimension(3,nz), intent(in) :: &
1949 : lhs_diff, & ! LHS diffustion terms
1950 : lhs_ma_zt ! LHS mean advection terms
1951 :
1952 : real( kind = core_rknd ), dimension(nz), intent(in) :: &
1953 : wind_speed, & ! wind speed; sqrt(u^2 + v^2) [m/s]
1954 : rho_ds_zm, & ! Dry, static density on momentum levels [kg/m^3]
1955 : invrs_rho_ds_zt ! Inv. dry, static density at thermo. levels [m^3/kg]
1956 :
1957 : real( kind = core_rknd ), intent(in) :: &
1958 : u_star_sqd ! Surface friction velocity, u_star, squared [m/s]
1959 :
1960 : logical, intent(in) :: &
1961 : l_imp_sfc_momentum_flux ! Flag for implicit momentum surface fluxes.
1962 :
1963 : type (stats_metadata_type), intent(in) :: &
1964 : stats_metadata
1965 :
1966 : !---------------------- InOut variables ----------------------
1967 : type (stats), target, intent(inout) :: &
1968 : stats_zt
1969 :
1970 : !---------------------- Local variables ----------------------
1971 : integer :: k, kp1, km1 ! Array indices
1972 :
1973 : real( kind = core_rknd ), dimension(nz) :: &
1974 0 : imp_sfc_flux
1975 :
1976 : ! Budget indices
1977 : integer :: ixm_ma, ixm_ta
1978 :
1979 : !---------------------- Begin Code ----------------------
1980 :
1981 0 : select case ( solve_type )
1982 : case ( windm_edsclrm_um )
1983 0 : ixm_ma = stats_metadata%ium_ma
1984 0 : ixm_ta = stats_metadata%ium_ta
1985 :
1986 : case ( windm_edsclrm_vm )
1987 0 : ixm_ma = stats_metadata%ivm_ma
1988 0 : ixm_ta = stats_metadata%ivm_ta
1989 :
1990 : case default
1991 0 : ixm_ma = 0
1992 0 : ixm_ta = 0
1993 :
1994 : end select
1995 :
1996 0 : imp_sfc_flux(:) = zero
1997 :
1998 0 : if ( l_imp_sfc_momentum_flux ) then
1999 :
2000 : ! Statistics: implicit contributions for um or vm.
2001 :
2002 : ! xm term ta is modified at level 2 to include the effects of the
2003 : ! surface flux. In this case, this effects the implicit portion of
2004 : ! the term, which handles the main diagonal for the turbulent advection term
2005 0 : if ( stats_metadata%ium_ta + stats_metadata%ivm_ta > 0 ) then
2006 0 : imp_sfc_flux(2) = - invrs_rho_ds_zt(2) * invrs_dzt(2) &
2007 0 : * rho_ds_zm(1) * ( u_star_sqd / wind_speed(2) )
2008 : endif
2009 :
2010 : endif ! l_imp_sfc_momentum_flux
2011 :
2012 : ! Finalize implicit contributions for xm
2013 :
2014 0 : do k = 2, nz-1, 1
2015 :
2016 0 : km1 = max( k-1, 1 )
2017 0 : kp1 = min( k+1, nz )
2018 :
2019 : ! xm mean advection
2020 : ! xm term ma is completely implicit; call stat_update_var_pt.
2021 : call stat_update_var_pt( ixm_ma, k, & ! intent(in)
2022 0 : - lhs_ma_zt(3,k) * xm(km1) &
2023 : - lhs_ma_zt(2,k) * xm(k) &
2024 0 : - lhs_ma_zt(1,k) * xm(kp1), & ! intent(in)
2025 0 : stats_zt ) ! intent(inout)
2026 :
2027 : ! xm turbulent transport (implicit component)
2028 : ! xm term ta has both implicit and explicit components;
2029 : ! call stat_end_update_pt.
2030 : call stat_end_update_pt( ixm_ta, k, & ! intent(in)
2031 0 : - 0.5_core_rknd * lhs_diff(3,k) * xm(km1) &
2032 : + ( - 0.5_core_rknd * lhs_diff(2,k) &
2033 : + imp_sfc_flux(k) ) * xm(k) &
2034 : - 0.5_core_rknd * lhs_diff(1,k) * xm(kp1), & ! intent(in)
2035 0 : stats_zt ) ! intent(inout)
2036 :
2037 : enddo
2038 :
2039 :
2040 : ! Upper boundary conditions
2041 0 : k = nz
2042 0 : km1 = max( k-1, 1 )
2043 :
2044 : ! xm mean advection
2045 : ! xm term ma is completely implicit; call stat_update_var_pt.
2046 : call stat_update_var_pt( ixm_ma, k, & ! intent(in)
2047 0 : - lhs_ma_zt(3,k) * xm(km1) &
2048 : - lhs_ma_zt(2,k) * xm(k), & ! intent(in)
2049 0 : stats_zt ) ! intent(inout)
2050 :
2051 : ! xm turbulent transport (implicit component)
2052 : ! xm term ta has both implicit and explicit components;
2053 : ! call stat_end_update_pt.
2054 : call stat_end_update_pt( ixm_ta, k, & ! intent(in)
2055 0 : - 0.5_core_rknd * lhs_diff(3,k) * xm(km1) &
2056 : + ( - 0.5_core_rknd * lhs_diff(2,k) &
2057 : + imp_sfc_flux(k) ) * xm(k), & ! intent(in)
2058 0 : stats_zt ) ! intent(inout)
2059 :
2060 :
2061 0 : return
2062 : end subroutine windm_edsclrm_implicit_stats
2063 :
2064 : !=============================================================================
2065 0 : subroutine compute_uv_tndcy( nz, ngrdcol, solve_type, &
2066 0 : fcor, perp_wind_m, perp_wind_g, &
2067 0 : xm_forcing, l_implemented, &
2068 : stats_metadata, &
2069 0 : stats_zt, &
2070 0 : xm_tndcy )
2071 :
2072 : ! Description:
2073 : ! Computes the explicit tendency for the um and vm wind components.
2074 : !
2075 : ! The only explicit tendency that is involved in the d(um)/dt or d(vm)/dt
2076 : ! equations is the Coriolis tendency.
2077 : !
2078 : ! The d(um)/dt equation contains the term:
2079 : !
2080 : ! - f * ( v_g - vm );
2081 : !
2082 : ! where f is the Coriolis parameter and v_g is the v component of the
2083 : ! geostrophic wind.
2084 : !
2085 : ! Likewise, the d(vm)/dt equation contains the term:
2086 : !
2087 : ! + f * ( u_g - um );
2088 : !
2089 : ! where u_g is the u component of the geostrophic wind.
2090 : !
2091 : ! This term is treated completely explicitly. The values of um, vm, u_g,
2092 : ! and v_g are all found on the thermodynamic levels.
2093 : !
2094 : ! Wind forcing from the GCSS cases is also added here.
2095 : !
2096 : ! References:
2097 : !-----------------------------------------------------------------------
2098 :
2099 : use grid_class, only: &
2100 : grid ! Type
2101 :
2102 : use stats_type_utilities, only: &
2103 : stat_update_var
2104 :
2105 : use stats_variables, only: &
2106 : stats_metadata_type
2107 :
2108 : use clubb_precision, only: &
2109 : core_rknd ! Variable(s)
2110 :
2111 : use stats_type, only: stats ! Type
2112 :
2113 : implicit none
2114 :
2115 : ! -------------------------- Input Variables --------------------------
2116 : integer, intent(in) :: &
2117 : nz, &
2118 : ngrdcol
2119 :
2120 : integer, intent(in) :: &
2121 : solve_type ! Description of what is being solved for
2122 :
2123 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
2124 : fcor ! Coriolis parameter [s^-1]
2125 :
2126 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
2127 : perp_wind_m, & ! Perpendicular component of the mean wind (e.g. v, for the u-eqn) [m/s]
2128 : perp_wind_g, & ! Perpendicular component of the geostropic wind (e.g. vg) [m/s]
2129 : xm_forcing ! Prescribed wind forcing [m/s/s]
2130 :
2131 : logical, intent(in) :: &
2132 : l_implemented ! Flag for CLUBB being implemented in a larger model.
2133 :
2134 : type (stats_metadata_type), intent(in) :: &
2135 : stats_metadata
2136 :
2137 : ! -------------------------- InOut Variables --------------------------
2138 : type (stats), target, dimension(ngrdcol), intent(inout) :: &
2139 : stats_zt
2140 :
2141 : ! -------------------------- Output Variables --------------------------
2142 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
2143 : xm_tndcy ! xm tendency [m/s^2]
2144 :
2145 : ! -------------------------- Local Variables --------------------------
2146 : integer :: &
2147 : ixm_gf, &
2148 : ixm_cf, &
2149 : ixm_f
2150 :
2151 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
2152 0 : xm_gf, &
2153 0 : xm_cf
2154 :
2155 : integer :: i, k
2156 :
2157 : ! -------------------------- Begin Code --------------------------
2158 :
2159 : !$acc enter data create( xm_gf, xm_cf )
2160 :
2161 0 : if ( .not. l_implemented ) then
2162 : ! Only compute the Coriolis term if the model is running on it's own,
2163 : ! and is not part of a larger, host model.
2164 :
2165 0 : select case ( solve_type )
2166 :
2167 : case ( windm_edsclrm_um )
2168 :
2169 0 : ixm_gf = stats_metadata%ium_gf
2170 0 : ixm_cf = stats_metadata%ium_cf
2171 0 : ixm_f = stats_metadata%ium_f
2172 :
2173 : !$acc parallel loop gang vector collapse(2) default(present)
2174 0 : do k = 1, nz
2175 0 : do i = 1, ngrdcol
2176 0 : xm_gf(i,k) = - fcor(i) * perp_wind_g(i,k)
2177 : end do
2178 : end do
2179 : !$acc end parallel loop
2180 :
2181 : !$acc parallel loop gang vector collapse(2) default(present)
2182 0 : do k = 1, nz
2183 0 : do i = 1, ngrdcol
2184 0 : xm_cf(i,k) = fcor(i) * perp_wind_m(i,k)
2185 : end do
2186 : end do
2187 : !$acc end parallel loop
2188 :
2189 : case ( windm_edsclrm_vm )
2190 :
2191 0 : ixm_gf = stats_metadata%ivm_gf
2192 0 : ixm_cf = stats_metadata%ivm_cf
2193 0 : ixm_f = stats_metadata%ivm_f
2194 :
2195 : !$acc parallel loop gang vector collapse(2) default(present)
2196 0 : do k = 1, nz
2197 0 : do i = 1, ngrdcol
2198 0 : xm_gf(i,k) = fcor(i) * perp_wind_g(i,k)
2199 : end do
2200 : end do
2201 : !$acc end parallel loop
2202 :
2203 : !$acc parallel loop gang vector collapse(2) default(present)
2204 0 : do k = 1, nz
2205 0 : do i = 1, ngrdcol
2206 0 : xm_cf(i,k) = -fcor(i) * perp_wind_m(i,k)
2207 : end do
2208 : end do
2209 : !$acc end parallel loop
2210 :
2211 : case default
2212 :
2213 0 : ixm_gf = 0
2214 0 : ixm_cf = 0
2215 0 : ixm_f = 0
2216 :
2217 : !$acc parallel loop gang vector collapse(2) default(present)
2218 0 : do k = 1, nz
2219 0 : do i = 1, ngrdcol
2220 0 : xm_gf(i,k) = 0._core_rknd
2221 0 : xm_cf(i,k) = 0._core_rknd
2222 : end do
2223 : end do
2224 : !$acc end parallel loop
2225 :
2226 : end select
2227 :
2228 : !$acc parallel loop gang vector collapse(2) default(present)
2229 0 : do k = 1, nz
2230 0 : do i = 1, ngrdcol
2231 0 : xm_tndcy(i,k) = xm_gf(i,k) + xm_cf(i,k) + xm_forcing(i,k)
2232 : end do
2233 : end do
2234 : !$acc end parallel loop
2235 :
2236 0 : if ( stats_metadata%l_stats_samp ) then
2237 :
2238 : !$acc update host( xm_gf, xm_cf, xm_forcing )
2239 :
2240 0 : do i = 1, ngrdcol
2241 : ! xm term gf is completely explicit; call stat_update_var.
2242 0 : call stat_update_var( ixm_gf, xm_gf(i,:), & ! intent(in)
2243 0 : stats_zt(i) ) ! intent(inout)
2244 :
2245 : ! xm term cf is completely explicit; call stat_update_var.
2246 : call stat_update_var( ixm_cf, xm_cf(i,:), & ! intent(in)
2247 0 : stats_zt(i) ) ! intent(inout)
2248 :
2249 : ! xm term F
2250 : call stat_update_var( ixm_f, xm_forcing(i,:), & ! intent(in)
2251 0 : stats_zt(i) ) ! intent(inout)
2252 : end do
2253 : endif
2254 :
2255 : else ! implemented in a host model.
2256 :
2257 : !$acc parallel loop gang vector collapse(2) default(present)
2258 0 : do k = 1, nz
2259 0 : do i = 1, ngrdcol
2260 0 : xm_tndcy(i,k) = 0.0_core_rknd
2261 : end do
2262 : end do
2263 : !$acc end parallel loop
2264 :
2265 : endif
2266 :
2267 : !$acc exit data delete( xm_gf, xm_cf )
2268 :
2269 0 : return
2270 :
2271 : end subroutine compute_uv_tndcy
2272 :
2273 : !======================================================================================
2274 352944 : subroutine windm_edsclrm_lhs( nz, ngrdcol, gr, dt, &
2275 352944 : lhs_ma_zt, lhs_diff, &
2276 352944 : wind_speed, u_star_sqd, &
2277 352944 : rho_ds_zm, invrs_rho_ds_zt, &
2278 : l_implemented, l_imp_sfc_momentum_flux, &
2279 352944 : lhs )
2280 : ! Description:
2281 : ! Calculate the implicit portion of the horizontal wind or eddy-scalar
2282 : ! time-tendency equation. See the description in subroutine
2283 : ! windm_edsclrm_solve for more details.
2284 : !
2285 : ! Notes:
2286 : ! Lower Boundary:
2287 : ! The lower boundary condition is a fixed-flux boundary condition, which
2288 : ! gets added into the time-tendency equation at level 2.
2289 : ! The value of xm(1) is located below the model surface and does not effect
2290 : ! the rest of the model. Since xm can be either a horizontal wind component
2291 : ! or a generic eddy scalar quantity, the value of xm(1) is simply set to the
2292 : ! value of xm(2) after the equation matrix has been solved.
2293 : !
2294 : ! --- THIS SUBROUTINE HAS BEEN OPTIMIZED ---
2295 : ! Simple changes to this procedure may adversely affect computational speed
2296 : ! - Gunther Huebler, Aug. 2018, clubb:ticket:834
2297 : !----------------------------------------------------------------------------------
2298 :
2299 : use grid_class, only: &
2300 : grid ! Type
2301 :
2302 : use clubb_precision, only: &
2303 : core_rknd ! Variable(s)
2304 :
2305 : implicit none
2306 :
2307 : ! ----------------------- Input Variables -----------------------
2308 : integer, intent(in) :: &
2309 : nz, &
2310 : ngrdcol
2311 :
2312 : type (grid), target, intent(in) :: gr
2313 :
2314 : real( kind = core_rknd ), intent(in) :: &
2315 : dt ! Model timestep [s]
2316 :
2317 : real( kind = core_rknd ), intent(in), dimension(3,ngrdcol,nz) :: &
2318 : lhs_diff, & ! LHS diffustion terms
2319 : lhs_ma_zt ! LHS mean advection terms
2320 :
2321 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
2322 : wind_speed, & ! wind speed; sqrt( u^2 + v^2 ) [m/s]
2323 : rho_ds_zm, & ! Dry, static density on momentum levels [kg/m^3]
2324 : invrs_rho_ds_zt ! Inv. dry, static density at thermo. levels [m^3/kg]
2325 :
2326 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
2327 : u_star_sqd ! Surface friction velocity, u_*, squared [m/s]
2328 :
2329 : logical, intent(in) :: &
2330 : l_implemented, & ! Flag for CLUBB being implemented in a larger model.
2331 : l_imp_sfc_momentum_flux ! Flag for implicit momentum surface fluxes.
2332 :
2333 : ! ----------------------- Output Variable -----------------------
2334 : real( kind = core_rknd ), dimension(3,ngrdcol,nz), intent(out) :: &
2335 : lhs ! Implicit contributions to xm (tridiagonal matrix)
2336 :
2337 : ! ----------------------- Local Variables -----------------------
2338 : real( kind = core_rknd ) :: &
2339 : invrs_dt ! Inverse of dt, 1/dt, used for computational efficiency
2340 :
2341 : integer :: k, i ! Loop variable
2342 :
2343 : ! ----------------------- Begin Code -----------------------
2344 :
2345 : ! Calculate coefs of eddy diffusivity and inverse of dt
2346 352944 : invrs_dt = 1.0_core_rknd / dt
2347 :
2348 : ! Set lower boundary, see notes
2349 : !$acc parallel loop gang vector default(present)
2350 5893344 : do i = 1, ngrdcol
2351 5540400 : lhs(1,i,1) = 0.0_core_rknd
2352 5540400 : lhs(2,i,1) = 1.0_core_rknd
2353 5893344 : lhs(3,i,1) = 0.0_core_rknd
2354 : end do
2355 : !$acc end parallel loop
2356 :
2357 : ! Add terms to lhs
2358 : !$acc parallel loop gang vector collapse(2) default(present)
2359 30000240 : do k = 2, nz
2360 495393840 : do i = 1, ngrdcol
2361 465393600 : lhs(1,i,k) = 0.5_core_rknd * lhs_diff(1,i,k)
2362 :
2363 465393600 : lhs(2,i,k) = 0.5_core_rknd * lhs_diff(2,i,k)
2364 :
2365 : ! LHS time tendency.
2366 465393600 : lhs(2,i,k) = lhs(2,i,k) + invrs_dt
2367 :
2368 495040896 : lhs(3,i,k) = 0.5_core_rknd * lhs_diff(3,i,k)
2369 : end do
2370 : end do
2371 : !$acc end parallel loop
2372 :
2373 : ! LHS mean advection term.
2374 352944 : if ( .not. l_implemented ) then
2375 : !$acc parallel loop gang vector collapse(2) default(present)
2376 0 : do k = 2, nz-1
2377 0 : do i = 1, ngrdcol
2378 0 : lhs(1:3,i,k) = lhs(1:3,i,k) + lhs_ma_zt(:,i,k)
2379 : end do
2380 : end do
2381 : !$acc end parallel loop
2382 :
2383 : endif
2384 :
2385 352944 : if ( l_imp_sfc_momentum_flux ) then
2386 :
2387 : ! LHS momentum surface flux.
2388 : !$acc parallel loop gang vector default(present)
2389 0 : do i = 1, ngrdcol
2390 0 : lhs(2,i,2) = lhs(2,i,2) + invrs_rho_ds_zt(i,2) * gr%invrs_dzt(i,2) &
2391 0 : * rho_ds_zm(i,1) * ( u_star_sqd(i) / wind_speed(i,2) )
2392 : end do
2393 : !$acc end parallel loop
2394 : end if ! l_imp_sfc_momentum_flux
2395 :
2396 352944 : return
2397 :
2398 : end subroutine windm_edsclrm_lhs
2399 :
2400 : !=============================================================================
2401 8117712 : subroutine windm_edsclrm_rhs( nz, ngrdcol, gr, solve_type, dt, &
2402 8117712 : lhs_diff, xm, xm_tndcy, &
2403 8117712 : rho_ds_zm, invrs_rho_ds_zt, &
2404 8117712 : l_imp_sfc_momentum_flux, xpwp_sfc, &
2405 : stats_metadata, &
2406 8117712 : stats_zt, &
2407 8117712 : rhs )
2408 : ! Description:
2409 : ! Calculate the explicit portion of the horizontal wind or eddy-scalar
2410 : ! time-tendency equation. See the description in subroutine
2411 : ! windm_edsclrm_solve for more details.
2412 : !
2413 : ! References:
2414 : ! None
2415 : !
2416 : ! Notes:
2417 : ! The lower boundary condition needs to be applied here at level 2.
2418 : ! The lower boundary condition is a "fixed flux" boundary condition.
2419 : ! The coding is the same as for a zero-flux boundary condition, but with
2420 : ! an extra term added on the right-hand side at the boundary level. For
2421 : ! the rest of the model code, a zero-flux boundary condition is applied
2422 : ! at level 1, and thus subroutine diffusion_zt_lhs is set-up to do that.
2423 : ! In order to apply the same boundary condition code here at level 2, an
2424 : ! adjuster needs to be used to tell diffusion_zt_lhs to use the code at
2425 : ! level 2 that it normally uses at level 1.
2426 : !
2427 : ! --- THIS SUBROUTINE HAS BEEN OPTIMIZED ---
2428 : ! Simple changes to this procedure may adversely affect computational speed
2429 : ! - Gunther Huebler, Aug. 2018, clubb:ticket:834
2430 : !----------------------------------------------------------------------------------
2431 :
2432 : use clubb_precision, only: &
2433 : core_rknd ! Variable(s)
2434 :
2435 : use diffusion, only: &
2436 : diffusion_zt_lhs ! Procedure(s)
2437 :
2438 : use stats_variables, only: &
2439 : stats_metadata_type
2440 :
2441 : use stats_type_utilities, only: &
2442 : stat_begin_update_pt, & ! Procedure(s)
2443 : stat_modify_pt
2444 :
2445 : use grid_class, only: &
2446 : grid ! Type
2447 :
2448 : use stats_type, only: stats ! Type
2449 :
2450 : implicit none
2451 :
2452 : !------------------- Input Variables -------------------
2453 : integer, intent(in) :: &
2454 : nz, &
2455 : ngrdcol
2456 :
2457 : type (grid), target, intent(in) :: &
2458 : gr
2459 :
2460 : integer, intent(in) :: &
2461 : solve_type ! Description of what is being solved for
2462 :
2463 : real( kind = core_rknd ), intent(in) :: &
2464 : dt ! Model timestep [s]
2465 :
2466 : real( kind = core_rknd ), dimension(3,ngrdcol,nz), intent(in) :: &
2467 : lhs_diff ! LHS diffustion terms
2468 :
2469 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
2470 : xm, & ! Eddy-scalar variable, xm (thermo. levels) [units vary]
2471 : xm_tndcy, & ! The explicit time-tendency acting on xm [units vary]
2472 : rho_ds_zm, & ! Dry, static density on momentum levels [kg/m^3]
2473 : invrs_rho_ds_zt ! Inv. dry, static density at thermo. levels [m^3/kg]
2474 :
2475 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
2476 : xpwp_sfc ! x'w' at the surface [units vary]
2477 :
2478 : logical, intent(in) :: &
2479 : l_imp_sfc_momentum_flux ! Flag for implicit momentum surface fluxes.
2480 :
2481 : type (stats_metadata_type), intent(in) :: &
2482 : stats_metadata
2483 :
2484 : !------------------- Inout Variable -------------------
2485 : type (stats), target, dimension(ngrdcol), intent(inout) :: &
2486 : stats_zt
2487 :
2488 : !------------------- Output Variable -------------------
2489 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
2490 : rhs ! Right-hand side (explicit) contributions.
2491 :
2492 : !------------------- Local Variables -------------------
2493 : integer :: i, k ! Loop variable
2494 :
2495 : real( kind = core_rknd ) :: invrs_dt
2496 :
2497 : integer :: ixm_ta
2498 :
2499 : !------------------- Begin Code -------------------
2500 :
2501 8117712 : invrs_dt = 1.0_core_rknd / dt ! Precalculate 1.0/dt to avoid redoing the divide
2502 :
2503 8117712 : select case ( solve_type )
2504 : case ( windm_edsclrm_um )
2505 0 : ixm_ta = stats_metadata%ium_ta
2506 : case ( windm_edsclrm_vm )
2507 0 : ixm_ta = stats_metadata%ivm_ta
2508 : case default ! Eddy scalars
2509 8117712 : ixm_ta = 0
2510 : end select
2511 :
2512 : ! For purposes of the matrix equation, rhs(1) is simply set to 0.
2513 : !$acc parallel loop gang vector default(present)
2514 135546912 : do i = 1, ngrdcol
2515 135546912 : rhs(i,1) = 0.0_core_rknd
2516 : end do
2517 : !$acc end parallel loop
2518 :
2519 : ! Non-boundary rhs calculation, this is a highly vectorized loop
2520 : !$acc parallel loop gang vector collapse(2) default(present)
2521 681887808 : do k = 2, nz-1
2522 11258511408 : do i = 1, ngrdcol
2523 21153247200 : rhs(i,k) = 0.5_core_rknd &
2524 10576623600 : * ( - lhs_diff(3,i,k) * xm(i,k-1) &
2525 : - lhs_diff(2,i,k) * xm(i,k) &
2526 10576623600 : - lhs_diff(1,i,k) * xm(i,k+1) ) &
2527 : + xm_tndcy(i,k) & ! RHS forcings
2528 53556888096 : + invrs_dt * xm(i,k) ! RHS time tendency
2529 : end do
2530 : end do
2531 : !$acc end parallel loop
2532 :
2533 : ! Upper boundary calculation
2534 : !$acc parallel loop gang vector default(present)
2535 135546912 : do i = 1, ngrdcol
2536 254858400 : rhs(i,nz) = 0.5_core_rknd &
2537 127429200 : * ( - lhs_diff(3,i,nz) * xm(i,nz-1) &
2538 : - lhs_diff(2,i,nz) * xm(i,nz) ) &
2539 : + xm_tndcy(i,nz) & ! RHS forcings
2540 390405312 : + invrs_dt * xm(i,nz) ! RHS time tendency
2541 : end do
2542 : !$acc end parallel loop
2543 :
2544 8117712 : if ( stats_metadata%l_stats_samp .and. ixm_ta > 0 ) then
2545 :
2546 : !$acc update host( lhs_diff, xm )
2547 :
2548 0 : do i = 1, ngrdcol
2549 :
2550 : ! Statistics: explicit contributions for um or vm.
2551 :
2552 : ! xm term ta has both implicit and explicit components; call
2553 : ! stat_begin_update_pt. Since stat_begin_update_pt automatically
2554 : ! subtracts the value sent in, reverse the sign on right-hand side
2555 : ! turbulent advection component.
2556 0 : do k = 2, nz-1
2557 :
2558 : call stat_begin_update_pt( ixm_ta, k, & ! intent(in)
2559 : 0.5_core_rknd &
2560 0 : * ( lhs_diff(3,i,k) * xm(i,k-1) &
2561 : + lhs_diff(2,i,k) * xm(i,k) & ! intent(in)
2562 0 : + lhs_diff(1,i,k) * xm(i,k+1) ), &
2563 0 : stats_zt(i) ) ! intent(inout)
2564 : end do
2565 :
2566 : ! Upper boundary
2567 : call stat_begin_update_pt( ixm_ta, nz, &
2568 : 0.5_core_rknd & ! intent(in)
2569 0 : * ( lhs_diff(3,i,nz) * xm(i,nz-1) &
2570 : + lhs_diff(2,i,nz) * xm(i,nz) ), & ! intent(in)
2571 0 : stats_zt(i) ) ! intent(inout)
2572 : end do
2573 : endif
2574 :
2575 8117712 : if ( .not. l_imp_sfc_momentum_flux ) then
2576 :
2577 : ! RHS generalized surface flux.
2578 : !$acc parallel loop gang vector default(present)
2579 135546912 : do i = 1, ngrdcol
2580 254858400 : rhs(i,2) = rhs(i,2) + invrs_rho_ds_zt(i,2) &
2581 0 : * gr%invrs_dzt(i,2) &
2582 262976112 : * rho_ds_zm(i,1) * xpwp_sfc(i)
2583 : end do
2584 : !$acc end parallel loop
2585 :
2586 8117712 : if ( stats_metadata%l_stats_samp .and. ixm_ta > 0 ) then
2587 :
2588 : !$acc update host( invrs_rho_ds_zt, rho_ds_zm, xpwp_sfc )
2589 :
2590 0 : do i = 1, ngrdcol
2591 : call stat_modify_pt( ixm_ta, 2, & ! intent(in)
2592 0 : + invrs_rho_ds_zt(i,2) &
2593 0 : * gr%invrs_dzt(i,2) &
2594 0 : * rho_ds_zm(i,1) * xpwp_sfc(i), & ! intent(in)
2595 0 : stats_zt(i) ) ! intent(inout)
2596 : end do
2597 : end if
2598 :
2599 : endif ! l_imp_sfc_momentum_flux
2600 :
2601 8117712 : return
2602 :
2603 : end subroutine windm_edsclrm_rhs
2604 :
2605 : end module advance_windm_edsclrm_module
|