Line data Source code
1 : !-----------------------------------------------------------------------
2 : ! $Id$
3 : !===============================================================================
4 : module diffusion
5 :
6 : ! Description:
7 : ! Module diffusion computes the eddy diffusion terms for all of the
8 : ! time-tendency (prognostic) equations in the CLUBB parameterization. Most of
9 : ! the eddy diffusion terms are solved for completely implicitly, and therefore
10 : ! become part of the left-hand side of their respective equations. However,
11 : ! wp2 and wp3 have an option to use a Crank-Nicholson eddy diffusion scheme,
12 : ! which has both implicit and explicit components.
13 : !
14 : ! Function diffusion_zt_lhs handles the eddy diffusion terms for the variables
15 : ! located at thermodynamic grid levels. These variables are: wp3 and all
16 : ! hydrometeor species. The variables um and vm also use the Crank-Nicholson
17 : ! eddy-diffusion scheme for their turbulent advection term.
18 : !
19 : ! Function diffusion_zm_lhs handles the eddy diffusion terms for the variables
20 : ! located at momentum grid levels. The variables are: wprtp, wpthlp, wp2,
21 : ! rtp2, thlp2, rtpthlp, up2, vp2, wpsclrp, sclrprtp, sclrpthlp, and sclrp2.
22 :
23 : implicit none
24 :
25 : private ! Default Scope
26 :
27 : public :: diffusion_zt_lhs, &
28 : diffusion_cloud_frac_zt_lhs, &
29 : diffusion_zm_lhs
30 :
31 : contains
32 :
33 : !=============================================================================
34 705888 : subroutine diffusion_zt_lhs( nz, ngrdcol, gr, K_zm, K_zt, nu, & ! In
35 705888 : invrs_rho_ds_zt, rho_ds_zm, & ! In
36 705888 : lhs ) ! Out
37 :
38 : ! Description:
39 : ! Vertical eddy diffusion of var_zt: implicit portion of the code.
40 : !
41 : ! The variable "var_zt" stands for a variable that is located at
42 : ! thermodynamic grid levels.
43 : !
44 : ! The d(var_zt)/dt equation contains an eddy diffusion term:
45 : !
46 : ! + d [ ( K_zm + nu ) * d(var_zt)/dz ] / dz.
47 : !
48 : ! This term is usually solved for completely implicitly, such that:
49 : !
50 : ! + d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz.
51 : !
52 : ! However, when a Crank-Nicholson scheme is used, the eddy diffusion term
53 : ! has both implicit and explicit components, such that:
54 : !
55 : ! + (1/2) * d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz
56 : ! + (1/2) * d [ ( K_zm + nu ) * d( var_zt(t) )/dz ] / dz;
57 : !
58 : ! for which the implicit component is:
59 : !
60 : ! + (1/2) * d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz.
61 : !
62 : ! Note: When the implicit term is brought over to the left-hand side,
63 : ! the sign is reversed and the leading "+" in front of the term
64 : ! is changed to a "-".
65 : !
66 : ! Timestep index (t) stands for the index of the current timestep, while
67 : ! timestep index (t+1) stands for the index of the next timestep, which is
68 : ! being advanced to in solving the d(var_zt)/dt equation.
69 : !
70 : ! The implicit portion of this term is discretized as follows:
71 : !
72 : ! The values of var_zt are found on the thermodynamic levels, while the
73 : ! values of K_zm are found on the momentum levels. The derivatives (d/dz)
74 : ! of var_zt are taken over the intermediate momentum levels. At the
75 : ! intermediate momentum levels, d(var_zt)/dz is multiplied by ( K_zm + nu ).
76 : ! Then, the derivative of the whole mathematical expression is taken over
77 : ! the central thermodynamic level, which yields the desired result.
78 : !
79 : ! --var_zt------------------------------------------------- t(k+1)
80 : !
81 : ! ==========d(var_zt)/dz==(K_zm+nu)======================== m(k)
82 : !
83 : ! --var_zt-------------------d[(K_zm+nu)*d(var_zt)/dz]/dz-- t(k)
84 : !
85 : ! ==========d(var_zt)/dz==(K_zm+nu)======================== m(k-1)
86 : !
87 : ! --var_zt------------------------------------------------- t(k-1)
88 : !
89 : ! The vertical indices t(k+1), m(k), t(k), m(k-1), and t(k-1) correspond
90 : ! with altitudes zt(k+1), zm(k), zt(k), zm(k-1), and zt(k-1), respectively.
91 : ! The letter "t" is used for thermodynamic levels and the letter "m" is used
92 : ! for momentum levels.
93 : !
94 : ! invrs_dzt(k) = 1 / ( zm(k) - zm(k-1) )
95 : ! invrs_dzm(k) = 1 / ( zt(k+1) - zt(k) )
96 : ! invrs_dzm(k-1) = 1 / ( zt(k) - zt(k-1) )
97 : !
98 : ! Note: This function only computes the general implicit form:
99 : ! + d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz.
100 : ! For a Crank-Nicholson scheme, the left-hand side result of this
101 : ! function will have to be multiplied by (1/2). For a
102 : ! Crank-Nicholson scheme, the right-hand side (explicit) component
103 : ! needs to be computed by multiplying the left-hand side results by
104 : ! (1/2), reversing the sign on each left-hand side element, and then
105 : ! multiplying each element by the appropriate var_zt(t) value from
106 : ! the appropriate vertical level.
107 : !
108 : !
109 : ! Boundary Conditions:
110 : !
111 : ! 1) Zero-flux boundary conditions.
112 : ! This function is set up to use zero-flux boundary conditions at both
113 : ! the lower boundary level and the upper boundary level. The flux, F,
114 : ! is the amount of var_zt flowing normal through the boundary per unit
115 : ! time per unit surface area. The derivative of the flux effects the
116 : ! time-tendency of var_zt, such that:
117 : !
118 : ! d(var_zt)/dt = -dF/dz.
119 : !
120 : ! For the 2nd-order eddy-diffusion term, +d[(K_zm+nu)*d(var_zt)/dz]/dz,
121 : ! the flux is:
122 : !
123 : ! F = -(K_zm+nu)*d(var_zt)/dz.
124 : !
125 : ! In order to have zero-flux boundary conditions, the derivative of
126 : ! var_zt, d(var_zt)/dz, needs to equal 0 at both the lower boundary and
127 : ! the upper boundary.
128 : !
129 : ! In order to discretize the lower boundary condition, consider a new
130 : ! level outside the model (thermodynamic level 0) just below the lower
131 : ! boundary level (thermodynamic level 1). The value of var_zt at the
132 : ! level just outside the model is defined to be the same as the value of
133 : ! var_zt at the lower boundary level. Therefore, the value of
134 : ! d(var_zt)/dz between the level just outside the model and the lower
135 : ! boundary level is 0, satisfying the zero-flux boundary condition. The
136 : ! other value for d(var_zt)/dz (between thermodynamic level 2 and
137 : ! thermodynamic level 1) is taken over the intermediate momentum level
138 : ! (momentum level 1), where it is multiplied by the factor
139 : ! ( K_zm(1) + nu ). Then, the derivative of the whole expression is
140 : ! taken over the central thermodynamic level.
141 : !
142 : ! -var_zt(2)-------------------------------------------- t(2)
143 : !
144 : ! ==========d(var_zt)/dz==(K_zm(1)+nu)================== m(1)
145 : !
146 : ! -var_zt(1)---------------d[(K_zm+nu)*d(var_zt)/dz]/dz- t(1) Boundary
147 : !
148 : ! [d(var_zt)/dz = 0]
149 : !
150 : ! -[var_zt(0) = var_zt(1)]-----(level outside model)---- t(0)
151 : !
152 : ! The result is dependent only on values of K_zm found at momentum
153 : ! level 1 and values of var_zt found at thermodynamic levels 1 and 2.
154 : ! Thus, it only affects 2 diagonals on the left-hand side matrix.
155 : !
156 : ! The same method can be used to discretize the upper boundary by
157 : ! considering a new level outside the model just above the upper boundary
158 : ! level.
159 : !
160 : ! 2) Fixed-point boundary conditions.
161 : ! Many equations in the model use fixed-point boundary conditions rather
162 : ! than zero-flux boundary conditions. This means that the value of
163 : ! var_zt stays the same over the course of the timestep at the lower
164 : ! boundary, as well as at the upper boundary.
165 : !
166 : ! In order to discretize the boundary conditions for equations requiring
167 : ! fixed-point boundary conditions, either:
168 : ! a) in the parent subroutine or function (that calls this function),
169 : ! loop over all vertical levels from the second-lowest to the
170 : ! second-highest, ignoring the boundary levels. Then set the values
171 : ! at the boundary levels in the parent subroutine; or
172 : ! b) in the parent subroutine or function, loop over all vertical levels
173 : ! and then overwrite the results at the boundary levels.
174 : !
175 : ! Either way, at the boundary levels, an array with a value of 1 at the
176 : ! main diagonal on the left-hand side and with values of 0 at all other
177 : ! diagonals on the left-hand side will preserve the right-hand side value
178 : ! at that level, thus satisfying the fixed-point boundary conditions.
179 : !
180 : !
181 : ! Conservation Properties:
182 : !
183 : ! When zero-flux boundary conditions are used, this technique of
184 : ! discretizing the eddy diffusion term leads to conservative differencing.
185 : ! When conservative differencing is in place, the column totals for each
186 : ! column in the left-hand side matrix (for the eddy diffusion term) should
187 : ! be equal to 0. This ensures that the total amount of the quantity var_zt
188 : ! over the entire vertical domain is being conserved, meaning that nothing
189 : ! is lost due to diffusional effects.
190 : !
191 : ! To see that this conservation law is satisfied, compute the eddy diffusion
192 : ! of var_zt and integrate vertically. In discretized matrix notation (where
193 : ! "i" stands for the matrix column and "j" stands for the matrix row):
194 : !
195 : ! 0 = Sum_j Sum_i ( 1/invrs_dzt )_i
196 : ! ( invrs_dzt * ((K_zm+nu)*invrs_dzm) )_ij (var_zt)_j.
197 : !
198 : ! The left-hand side matrix, ( invrs_dzt * ((K_zm+nu)*invrs_dzm) )_ij, is
199 : ! partially written below. The sum over i in the above equation removes
200 : ! invrs_dzt everywhere from the matrix below. The sum over j leaves the
201 : ! column totals that are desired.
202 : !
203 : ! Left-hand side matrix contributions from eddy diffusion term; first four
204 : ! vertical levels:
205 : !
206 : ! -------------------------------------------------------------------------->
207 : !k=1 | +invrs_dzt(k) -invrs_dzt(k) 0
208 : ! | *(K_zm(k)+nu) *(K_zm(k)+nu)
209 : ! | *invrs_dzm(k) *invrs_dzm(k)
210 : ! |
211 : !k=2 | -invrs_dzt(k) +invrs_dzt(k) -invrs_dzt(k)
212 : ! | *(K_zm(k-1)+nu) *[ (K_zm(k)+nu) *(K_zm(k)+nu)
213 : ! | *invrs_dzm(k-1) *invrs_dzm(k) *invrs_dzm(k)
214 : ! | +(K_zm(k-1)+nu)
215 : ! | *invrs_dzm(k-1) ]
216 : ! |
217 : !k=3 | 0 -invrs_dzt(k) +invrs_dzt(k)
218 : ! | *(K_zm(k-1)+nu) *[ (K_zm(k)+nu)
219 : ! | *invrs_dzm(k-1) *invrs_dzm(k)
220 : ! | +(K_zm(k-1)+nu)
221 : ! | *invrs_dzm(k-1) ]
222 : ! |
223 : !k=4 | 0 0 -invrs_dzt(k)
224 : ! | *(K_zm(k-1)+nu)
225 : ! | *invrs_dzm(k-1)
226 : ! \ /
227 : !
228 : ! Note: The superdiagonal term from level 3 and both the main diagonal and
229 : ! superdiagonal terms from level 4 are not shown on this diagram.
230 : !
231 : ! Note: The matrix shown is a tridiagonal matrix. For a band diagonal
232 : ! matrix (with 5 diagonals), there would be an extra row between each
233 : ! of the rows shown and an extra column between each of the columns
234 : ! shown. However, for the purposes of the var_zt eddy diffusion
235 : ! term, those extra row and column values are all 0, and the
236 : ! conservation properties of the matrix aren't effected.
237 : !
238 : ! If fixed-point boundary conditions are used, the matrix entries at
239 : ! level 1 (k=1) read: 1 0 0; which means that conservative differencing
240 : ! is not in play. The total amount of var_zt over the entire vertical
241 : ! domain is not being conserved, as amounts of var_zt may be fluxed out
242 : ! through the upper boundary or lower boundary through the effects of
243 : ! diffusion.
244 : !
245 : ! Brian Griffin. April 26, 2008.
246 :
247 : ! References:
248 : ! None
249 : !-----------------------------------------------------------------------
250 :
251 : use grid_class, only: &
252 : grid, & ! Type
253 : ddzm ! Procedure
254 :
255 : use constants_clubb, only: &
256 : zero ! Constant(s)
257 :
258 : use clubb_precision, only: &
259 : core_rknd ! Variable(s)
260 :
261 : use model_flags, only: &
262 : l_upwind_Kh_dp_term
263 :
264 : implicit none
265 :
266 : ! Constant parameters
267 : integer, parameter :: &
268 : kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
269 : k_tdiag = 2, & ! Thermodynamic main diagonal index.
270 : km1_tdiag = 3 ! Thermodynamic subdiagonal index.
271 :
272 : !------------------------ Input Variables ------------------------
273 : integer, intent(in) :: &
274 : nz, &
275 : ngrdcol
276 :
277 : type (grid), target, intent(in) :: gr
278 :
279 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
280 : K_zm, & ! Coef. of eddy diffusivity at momentum levels [m^2/s]
281 : K_zt, & ! Coef. of eddy diffusivity at thermo. levels [m^2/s]
282 : rho_ds_zm, & ! Dry statis density on momentum levels [kg]
283 : invrs_rho_ds_zt ! Inverse dry statis density on thermo. levels [1/kg]
284 :
285 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
286 : nu ! Background constant coef. of eddy diffusivity [m^2/s]
287 :
288 : !------------------------ Return Variable ------------------------
289 : real( kind = core_rknd ), dimension(3,ngrdcol,nz), intent(out) :: &
290 : lhs ! LHS coefficient of diffusion term [1/s]
291 :
292 : !------------------------ Local Variable ------------------------
293 : integer :: i,k ! Vertical level index
294 :
295 : real( kind = core_rknd ), dimension(3,ngrdcol,nz) :: &
296 705888 : lhs_upwind ! LHS coefficient diffusion term due to upwind [1/s]
297 :
298 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
299 705888 : drhoKdz_zt, &
300 1411776 : K_zm_nu, &
301 705888 : rho_K_zm_nu, &
302 705888 : ddzm_rho_K_zm_nu
303 :
304 : !------------------------ Begin code ------------------------
305 :
306 : !$acc data copyin( gr, gr%invrs_dzm, gr%invrs_dzt, &
307 : !$acc K_zm, K_zt, rho_ds_zm, invrs_rho_ds_zt, nu ) &
308 : !$acc copyout( lhs ) &
309 : !$acc create( lhs_upwind, drhoKdz_zt, K_zm_nu, rho_K_zm_nu, ddzm_rho_K_zm_nu )
310 :
311 : !$acc parallel loop gang vector collapse(2) default(present)
312 60706368 : do k = 1, nz
313 1002574368 : do i = 1, ngrdcol
314 1001868480 : K_zm_nu(i,k) = K_zm(i,k) + nu(i)
315 : end do
316 : end do
317 : !$acc end parallel loop
318 :
319 : if ( l_upwind_Kh_dp_term ) then
320 :
321 : !$acc parallel loop gang vector collapse(2) default(present)
322 : do k = 1, nz
323 : do i = 1, ngrdcol
324 : rho_K_zm_nu(i,k) = rho_ds_zm(i,k) * K_zm_nu(i,k)
325 : end do
326 : end do
327 : !$acc end parallel loop
328 :
329 : ! calculate the dKh_zt/dz
330 : ddzm_rho_K_zm_nu = ddzm( nz, ngrdcol, gr, rho_K_zm_nu )
331 :
332 : !$acc parallel loop gang vector collapse(2) default(present)
333 : do k = 1, nz
334 : do i = 1, ngrdcol
335 : drhoKdz_zt(i,k) = - invrs_rho_ds_zt(i,k) * ddzm_rho_K_zm_nu(i,k)
336 : end do
337 : end do
338 : !$acc end parallel loop
339 :
340 :
341 : ! extra terms with upwind scheme
342 : ! k = 1 (bottom level); lowere boundary level
343 : !$acc parallel loop gang vector default(present)
344 : do i = 1, ngrdcol
345 : lhs_upwind(kp1_tdiag,i,1) = + min( drhoKdz_zt(i,1) , zero ) * gr%invrs_dzm(i,1)
346 : lhs_upwind(k_tdiag,i,1) = - min( drhoKdz_zt(i,1) , zero ) * gr%invrs_dzm(i,1)
347 : lhs_upwind(km1_tdiag,i,1) = zero
348 : end do
349 : !$acc end parallel loop
350 :
351 : ! Most of the interior model; normal conditions.
352 : !$acc parallel loop gang vector collapse(2) default(present)
353 : do k = 2, nz-1, 1
354 : do i = 1, ngrdcol
355 :
356 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
357 : lhs_upwind(kp1_tdiag,i,k) = + min( drhoKdz_zt(i,k) , zero ) * gr%invrs_dzm(i,k)
358 :
359 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
360 : lhs_upwind(k_tdiag,i,k) = - min( drhoKdz_zt(i,k) , zero ) * gr%invrs_dzm(i,k) &
361 : + max( drhoKdz_zt(i,k) , zero ) * gr%invrs_dzm(i,k-1)
362 :
363 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
364 : lhs_upwind(km1_tdiag,i,k) = - max( drhoKdz_zt(i,k) , zero ) * gr%invrs_dzm(i,k-1)
365 :
366 : end do
367 : end do
368 : !$acc end parallel loop
369 :
370 : ! k = nz (top level); upper boundary level.
371 : ! Only relevant if zero-flux boundary conditions are used.
372 : !$acc parallel loop gang vector default(present)
373 : do i = 1, ngrdcol
374 : lhs_upwind(kp1_tdiag,i,nz) = zero
375 : lhs_upwind(k_tdiag,i,nz) = + max( drhoKdz_zt(i,nz) , zero ) * gr%invrs_dzm(i,nz-1)
376 : lhs_upwind(km1_tdiag,i,nz) = - max( drhoKdz_zt(i,nz) , zero ) * gr%invrs_dzm(i,nz-1)
377 : end do
378 : !$acc end parallel loop
379 :
380 : end if
381 :
382 : ! k = 1 (bottom level); lower boundary level.
383 : ! Only relevant if zero-flux boundary conditions are used.
384 : ! These k=1 lines currently do not have any effect on model results.
385 : ! These k=1 level of this "lhs" array is not fed into
386 : ! the final LHS matrix that will be used to solve for the next timestep.
387 :
388 : if ( .not. l_upwind_Kh_dp_term ) then
389 :
390 : !$acc parallel loop gang vector default(present)
391 11786688 : do i = 1, ngrdcol
392 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
393 22161600 : lhs(kp1_tdiag,i,1) = - gr%invrs_dzt(i,1) * invrs_rho_ds_zt(i,1) &
394 22161600 : * ( K_zm(i,1) + nu(i) ) * rho_ds_zm(i,1) * gr%invrs_dzm(i,1)
395 :
396 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
397 0 : lhs(k_tdiag,i,1) = + gr%invrs_dzt(i,1) * invrs_rho_ds_zt(i,1) &
398 11080800 : * ( K_zm(i,1) + nu(i) ) * rho_ds_zm(i,1) * gr%invrs_dzm(i,1)
399 :
400 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
401 11786688 : lhs(km1_tdiag,i,1) = zero
402 : end do
403 : !$acc end parallel loop
404 :
405 : else
406 :
407 : !$acc parallel loop gang vector default(present)
408 : do i = 1, ngrdcol
409 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
410 : lhs(kp1_tdiag,i,1) = - gr%invrs_dzt(i,1) * ( K_zt(i,1) + nu(i) ) * gr%invrs_dzm(i,1) &
411 : + lhs_upwind(kp1_tdiag,i,1)
412 :
413 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
414 : lhs(k_tdiag,i,1) = + gr%invrs_dzt(i,1) * ( K_zt(i,1) + nu(i) ) * gr%invrs_dzm(i,1) &
415 : + lhs_upwind(k_tdiag,i,1)
416 :
417 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
418 : lhs(km1_tdiag,i,1) = zero + lhs_upwind(km1_tdiag,i,1)
419 : end do
420 : !$acc end parallel loop
421 :
422 : end if
423 :
424 : if ( .not. l_upwind_Kh_dp_term ) then
425 :
426 : ! Most of the interior model; normal conditions.
427 : !$acc parallel loop gang vector collapse(2) default(present)
428 59294592 : do k = 2, nz-1, 1
429 979000992 : do i = 1, ngrdcol
430 :
431 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
432 1839412800 : lhs(kp1_tdiag,i,k) &
433 0 : = - gr%invrs_dzt(i,k) * invrs_rho_ds_zt(i,k) &
434 2759119200 : * K_zm_nu(i,k) * rho_ds_zm(i,k) * gr%invrs_dzm(i,k)
435 :
436 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
437 : lhs(k_tdiag,i,k) &
438 0 : = + gr%invrs_dzt(i,k) * invrs_rho_ds_zt(i,k) &
439 0 : * ( K_zm_nu(i,k) * rho_ds_zm(i,k) * gr%invrs_dzm(i,k) &
440 919706400 : + K_zm_nu(i,k-1) * rho_ds_zm(i,k-1) * gr%invrs_dzm(i,k-1) )
441 :
442 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
443 : lhs(km1_tdiag,i,k) &
444 0 : = - gr%invrs_dzt(i,k) * invrs_rho_ds_zt(i,k) &
445 978295104 : * K_zm_nu(i,k-1) * rho_ds_zm(i,k-1) * gr%invrs_dzm(i,k-1)
446 : end do
447 : end do ! k = 2, nz-1, 1
448 : !$acc end parallel loop
449 :
450 : else
451 :
452 : !$acc parallel loop gang vector collapse(2) default(present)
453 : do k = 2, nz-1, 1
454 : do i = 1, ngrdcol
455 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
456 : lhs(kp1_tdiag,i,k) &
457 : = - gr%invrs_dzt(i,k) * ( K_zt(i,k) + nu(i) ) * gr%invrs_dzm(i,k) &
458 : + lhs_upwind(kp1_tdiag,i,k)
459 :
460 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
461 : lhs(k_tdiag,i,k) &
462 : = + gr%invrs_dzt(i,k) * ( K_zt(i,k) + nu(i) ) * ( gr%invrs_dzm(i,k) + gr%invrs_dzm(i,k-1) ) &
463 : + lhs_upwind(k_tdiag,i,k)
464 :
465 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
466 : lhs(km1_tdiag,i,k) &
467 : = - gr%invrs_dzt(i,k) * ( K_zt(i,k) + nu(i) ) * gr%invrs_dzm(i,k-1) &
468 : + lhs_upwind(km1_tdiag,i,k)
469 : end do
470 : end do ! k = 2, nz-1, 1
471 : !$acc end parallel loop
472 :
473 : end if
474 :
475 : ! k = nz (top level); upper boundary level.
476 : ! Only relevant if zero-flux boundary conditions are used.
477 :
478 : if ( .not. l_upwind_Kh_dp_term ) then
479 :
480 : !$acc parallel loop gang vector default(present)
481 11786688 : do i = 1, ngrdcol
482 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
483 11080800 : lhs(kp1_tdiag,i,nz) = zero
484 :
485 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
486 : lhs(k_tdiag,i,nz) &
487 0 : = + gr%invrs_dzt(i,nz) * invrs_rho_ds_zt(i,nz) &
488 11080800 : * ( K_zm(i,nz-1) + nu(i) ) * rho_ds_zm(i,nz-1) * gr%invrs_dzm(i,nz-1)
489 :
490 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
491 : lhs(km1_tdiag,i,nz) &
492 0 : = - gr%invrs_dzt(i,nz) * invrs_rho_ds_zt(i,nz) &
493 11786688 : * ( K_zm(i,nz-1) + nu(i) ) * rho_ds_zm(i,nz-1) * gr%invrs_dzm(i,nz-1)
494 : end do
495 : !$acc end parallel loop
496 :
497 : else
498 :
499 : !$acc parallel loop gang vector default(present)
500 : do i = 1, ngrdcol
501 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
502 : lhs(kp1_tdiag,i,nz) = zero + lhs_upwind(kp1_tdiag,i,nz)
503 :
504 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
505 : lhs(k_tdiag,i,nz) &
506 : = + gr%invrs_dzt(i,nz) * ( K_zt(i,nz) + nu(i) ) * gr%invrs_dzm(i,nz-1) &
507 : + lhs_upwind(k_tdiag,i,nz)
508 :
509 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
510 : lhs(km1_tdiag,i,nz) &
511 : = - gr%invrs_dzt(i,nz) * ( K_zt(i,nz) + nu(i) ) * gr%invrs_dzm(i,nz-1) &
512 : + lhs_upwind(km1_tdiag,i,nz)
513 : end do
514 : !$acc end parallel loop
515 :
516 : end if
517 :
518 : !$acc end data
519 :
520 705888 : return
521 :
522 : end subroutine diffusion_zt_lhs
523 :
524 : !=============================================================================
525 0 : function diffusion_cloud_frac_zt_lhs &
526 : ( gr, K_zm, K_zmm1, cloud_frac_zt, cloud_frac_ztm1, &
527 : cloud_frac_ztp1, cloud_frac_zm, &
528 : cloud_frac_zmm1, &
529 : nu, invrs_dzmm1, invrs_dzm, invrs_dzt, level ) &
530 0 : result( lhs )
531 :
532 : ! Description:
533 : ! This function adds a weight of cloud fraction to the existing diffusion
534 : ! function for number concentration variables (e.g. cloud droplet number
535 : ! concentration). This code should be considered experimental and may
536 : ! contain bugs.
537 : ! References:
538 : ! This algorithm uses equations derived from Guo, et al. 2009.
539 : !-----------------------------------------------------------------------------
540 :
541 : use grid_class, only: &
542 : grid ! Type
543 :
544 : use clubb_precision, only: &
545 : core_rknd ! Variable(s)
546 :
547 : implicit none
548 :
549 : type (grid), target, intent(in) :: gr
550 :
551 : ! External
552 : intrinsic :: min
553 :
554 : ! Constant parameters
555 : real( kind = core_rknd ), parameter :: &
556 : cf_ratio = 10._core_rknd ! Maximum cloud-fraction coefficient applied to Kh_zm
557 :
558 : integer, parameter :: &
559 : kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
560 : k_tdiag = 2, & ! Thermodynamic main diagonal index.
561 : km1_tdiag = 3 ! Thermodynamic subdiagonal index.
562 :
563 : ! Input Variables
564 : real( kind = core_rknd ), intent(in) :: &
565 : K_zm, & ! Coef. of eddy diffusivity at mom. level (k) [m^2/s]
566 : K_zmm1, & ! Coef. of eddy diffusivity at mom. level (k-1) [m^2/s]
567 : cloud_frac_zt, & ! Cloud fraction at the thermo. level (k) [-]
568 : cloud_frac_ztm1, & ! Cloud fraction at the thermo. level (k-1) [-]
569 : cloud_frac_ztp1, & ! Cloud fraction at the thermo. level (k+1) [-]
570 : cloud_frac_zm, & ! Cloud fraction at the momentum level (k) [-]
571 : cloud_frac_zmm1, & ! Cloud fraction at the momentum level (k-1) [-]
572 : invrs_dzt, & ! Inverse of grid spacing over thermo. lev. (k) [1/m]
573 : invrs_dzm, & ! Inverse of grid spacing over mom. level (k) [1/m]
574 : invrs_dzmm1 ! Inverse of grid spacing over mom. level (k-1) [1/m]
575 :
576 : real( kind = core_rknd ), intent(in) :: &
577 : nu ! Background constant coef. of eddy diffusivity [m^2/s]
578 :
579 : integer, intent(in) :: &
580 : level ! Thermodynamic level where calculation occurs. [-]
581 :
582 : ! Return Variable
583 : real( kind = core_rknd ), dimension(3) :: lhs
584 :
585 : ! ---- Begin Code ----
586 :
587 0 : if ( level == 1 ) then
588 :
589 : ! k = 1 (bottom level); lower boundary level.
590 : ! Only relevant if zero-flux boundary conditions are used.
591 :
592 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
593 : ! lhs(kp1_tdiag) = - invrs_dzt &
594 : ! * (K_zm+nu) &
595 : ! * ( cloud_frac_zm / cloud_frac_ztp1 ) * invrs_dzm
596 : lhs(kp1_tdiag) = - invrs_dzt &
597 : * (K_zm &
598 : * min( cloud_frac_zm / cloud_frac_ztp1, cf_ratio ) &
599 0 : + nu) * invrs_dzm
600 :
601 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
602 : ! lhs(k_tdiag) = + invrs_dzt &
603 : ! * (K_zm+nu) &
604 : ! * ( cloud_frac_zm / cloud_frac_ztp1 ) * invrs_dzm
605 : lhs(k_tdiag) = + invrs_dzt &
606 : * (K_zm &
607 : * min( cloud_frac_zm / cloud_frac_ztp1, cf_ratio ) &
608 0 : + nu) * invrs_dzm
609 :
610 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
611 0 : lhs(km1_tdiag) = 0.0_core_rknd
612 :
613 :
614 0 : else if ( level > 1 .and. level < gr%nz ) then
615 :
616 : ! Most of the interior model; normal conditions.
617 :
618 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
619 : ! lhs(kp1_tdiag) = - invrs_dzt &
620 : ! * (K_zm+nu) &
621 : ! * ( cloud_frac_zm / cloud_frac_ztp1 ) * invrs_dzm
622 : ! lhs(kp1_tdiag) = - invrs_dzt &
623 : ! * (K_zm &
624 : ! * ( cloud_frac_zm / cloud_frac_ztp1 ) &
625 : ! + nu ) * invrs_dzm
626 : lhs(kp1_tdiag) = - invrs_dzt &
627 : * (K_zm &
628 : * min( cloud_frac_zm / cloud_frac_ztp1, cf_ratio ) &
629 0 : + nu ) * invrs_dzm
630 :
631 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
632 : ! lhs(k_tdiag) = + invrs_dzt &
633 : ! * ( ((K_zm+nu)*cloud_frac_zm)*invrs_dzm &
634 : ! + ((K_zmm1+nu)*cloud_frac_zmm1)*invrs_dzmm1 ) &
635 : ! / cloud_frac_zt
636 : ! lhs(k_tdiag) = + invrs_dzt &
637 : ! * ( nu*(invrs_dzm+invrs_dzmm1) + &
638 : ! ( ((K_zm*cloud_frac_zm)*invrs_dzm +
639 : ! (K_zmm1*cloud_frac_zmm1)*invrs_dzmm1)&
640 : ! / cloud_frac_zt &
641 : ! ) &
642 : ! )
643 : lhs(k_tdiag) = + invrs_dzt &
644 : * ( nu*(invrs_dzm+invrs_dzmm1) + &
645 : ( K_zm*invrs_dzm* &
646 : min( cloud_frac_zm / cloud_frac_zt, &
647 : cf_ratio ) &
648 : + K_zmm1*invrs_dzmm1* &
649 : min( cloud_frac_zmm1 / cloud_frac_zt, &
650 : cf_ratio ) &
651 : ) &
652 0 : )
653 :
654 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
655 : ! lhs(km1_tdiag) = - invrs_dzt * (K_zmm1+nu) * &
656 : ! ( cloud_frac_zmm1 / cloud_frac_ztm1 ) * invrs_dzmm1
657 : lhs(km1_tdiag) = - invrs_dzt &
658 : * (K_zmm1 &
659 : * min( cloud_frac_zmm1 / cloud_frac_ztm1, &
660 : cf_ratio ) &
661 0 : + nu ) * invrs_dzmm1
662 :
663 0 : else if ( level == gr%nz ) then
664 :
665 : ! k = gr%nz (top level); upper boundary level.
666 : ! Only relevant if zero-flux boundary conditions are used.
667 :
668 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
669 0 : lhs(kp1_tdiag) = 0.0_core_rknd
670 :
671 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
672 : ! lhs(k_tdiag) = + invrs_dzt &
673 : ! *(K_zmm1+nu) &
674 : ! *( cloud_frac_zmm1 / cloud_frac_ztm1 ) * invrs_dzmm1
675 : lhs(k_tdiag) = + invrs_dzt &
676 : * (K_zmm1 &
677 : * min( cloud_frac_zmm1 / cloud_frac_ztm1, &
678 : cf_ratio ) &
679 0 : + nu) * invrs_dzmm1
680 :
681 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
682 : ! lhs(km1_tdiag) = - invrs_dzt * (K_zmm1+nu) * &
683 : ! ( cloud_frac_zmm1 / cloud_frac_ztm1 ) * invrs_dzmm1
684 : lhs(km1_tdiag) = - invrs_dzt &
685 : * (K_zmm1 &
686 : * min( cloud_frac_zmm1 / cloud_frac_ztm1, &
687 : cf_ratio ) &
688 0 : + nu) * invrs_dzmm1
689 :
690 : end if
691 :
692 0 : return
693 0 : end function diffusion_cloud_frac_zt_lhs
694 :
695 : !=============================================================================
696 1411776 : subroutine diffusion_zm_lhs( nz, ngrdcol, gr, K_zt, K_zm, nu, & ! In
697 1411776 : invrs_rho_ds_zm, rho_ds_zt, & ! In
698 1411776 : lhs ) ! Out
699 :
700 : ! Description:
701 : ! Vertical eddy diffusion of var_zm: implicit portion of the code.
702 : !
703 : ! The variable "var_zm" stands for a variable that is located at momentum
704 : ! grid levels.
705 : !
706 : ! The d(var_zm)/dt equation contains an eddy diffusion term:
707 : !
708 : ! + d [ ( K_zt + nu ) * d(var_zm)/dz ] / dz.
709 : !
710 : ! This term is usually solved for completely implicitly, such that:
711 : !
712 : ! + d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz.
713 : !
714 : ! However, when a Crank-Nicholson scheme is used, the eddy diffusion term
715 : ! has both implicit and explicit components, such that:
716 : !
717 : ! + (1/2) * d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz
718 : ! + (1/2) * d [ ( K_zt + nu ) * d( var_zm(t) )/dz ] / dz;
719 : !
720 : ! for which the implicit component is:
721 : !
722 : ! + (1/2) * d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz.
723 : !
724 : ! Note: When the implicit term is brought over to the left-hand side,
725 : ! the sign is reversed and the leading "+" in front of the term
726 : ! is changed to a "-".
727 : !
728 : ! Timestep index (t) stands for the index of the current timestep, while
729 : ! timestep index (t+1) stands for the index of the next timestep, which is
730 : ! being advanced to in solving the d(var_zm)/dt equation.
731 : !
732 : ! The implicit portion of this term is discretized as follows:
733 : !
734 : ! The values of var_zm are found on the momentum levels, while the values of
735 : ! K_zt are found on the thermodynamic levels. The derivatives (d/dz) of
736 : ! var_zm are taken over the intermediate thermodynamic levels. At the
737 : ! intermediate thermodynamic levels, d(var_zm)/dz is multiplied by
738 : ! ( K_zt + nu ). Then, the derivative of the whole mathematical expression
739 : ! is taken over the central momentum level, which yields the desired result.
740 : !
741 : ! ==var_zm================================================= m(k+1)
742 : !
743 : ! ----------d(var_zm)/dz--(K_zt+nu)------------------------ t(k+1)
744 : !
745 : ! ==var_zm===================d[(K_zt+nu)*d(var_zm)/dz]/dz== m(k)
746 : !
747 : ! ----------d(var_zm)/dz--(K_zt+nu)------------------------ t(k)
748 : !
749 : ! ==var_zm================================================= m(k-1)
750 : !
751 : ! The vertical indices m(k+1), t(k+1), m(k), t(k), and m(k-1) correspond
752 : ! with altitudes zm(k+1), zt(k+1), zm(k), zt(k), and zm(k-1), respectively.
753 : ! The letter "t" is used for thermodynamic levels and the letter "m" is used
754 : ! for momentum levels.
755 : !
756 : ! invrs_dzm(k) = 1 / ( zt(k+1) - zt(k) )
757 : ! invrs_dzt(k+1) = 1 / ( zm(k+1) - zm(k) )
758 : ! invrs_dzt(k) = 1 / ( zm(k) - zm(k-1) )
759 : !
760 : ! Note: This function only computes the general implicit form:
761 : ! + d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz.
762 : ! For a Crank-Nicholson scheme, the left-hand side result of this
763 : ! function will have to be multiplied by (1/2). For a
764 : ! Crank-Nicholson scheme, the right-hand side (explicit) component
765 : ! needs to be computed by multiplying the left-hand side results by
766 : ! (1/2), reversing the sign on each left-hand side element, and then
767 : ! multiplying each element by the appropriate var_zm(t) value from
768 : ! the appropriate vertical level.
769 : !
770 : !
771 : ! Boundary Conditions:
772 : !
773 : ! 1) Zero-flux boundary conditions.
774 : ! This function is set up to use zero-flux boundary conditions at both
775 : ! the lower boundary level and the upper boundary level. The flux, F,
776 : ! is the amount of var_zm flowing normal through the boundary per unit
777 : ! time per unit surface area. The derivative of the flux effects the
778 : ! time-tendency of var_zm, such that:
779 : !
780 : ! d(var_zm)/dt = -dF/dz.
781 : !
782 : ! For the 2nd-order eddy-diffusion term, +d[(K_zt+nu)*d(var_zm)/dz]/dz,
783 : ! the flux is:
784 : !
785 : ! F = -(K_zt+nu)*d(var_zm)/dz.
786 : !
787 : ! In order to have zero-flux boundary conditions, the derivative of
788 : ! var_zm, d(var_zm)/dz, needs to equal 0 at both the lower boundary and
789 : ! the upper boundary.
790 : !
791 : ! In order to discretize the lower boundary condition, consider a new
792 : ! level outside the model (momentum level 0) just below the lower
793 : ! boundary level (momentum level 1). The value of var_zm at the level
794 : ! just outside the model is defined to be the same as the value of var_zm
795 : ! at the lower boundary level. Therefore, the value of d(var_zm)/dz
796 : ! between the level just outside the model and the lower boundary level
797 : ! is 0, satisfying the zero-flux boundary condition. The other value for
798 : ! d(var_zm)/dz (between momentum level 2 and momentum level 1) is taken
799 : ! over the intermediate thermodynamic level (thermodynamic level 2),
800 : ! where it is multiplied by the factor ( K_zt(2) + nu ). Then, the
801 : ! derivative of the whole expression is taken over the central momentum
802 : ! level.
803 : !
804 : ! =var_zm(2)============================================ m(2)
805 : !
806 : ! ----------d(var_zm)/dz==(K_zt(2)+nu)------------------ t(2)
807 : !
808 : ! =var_zm(1)===============d[(K_zt+nu)*d(var_zm)/dz]/dz= m(1) Boundary
809 : !
810 : ! ----------[d(var_zm)/dz = 0]-------------------------- t(1)
811 : !
812 : ! =[var_zm(0) = var_zm(1)]=====(level outside model)==== m(0)
813 : !
814 : ! The result is dependent only on values of K_zt found at thermodynamic
815 : ! level 2 and values of var_zm found at momentum levels 1 and 2. Thus,
816 : ! it only affects 2 diagonals on the left-hand side matrix.
817 : !
818 : ! The same method can be used to discretize the upper boundary by
819 : ! considering a new level outside the model just above the upper boundary
820 : ! level.
821 : !
822 : ! 2) Fixed-point boundary conditions.
823 : ! Many equations in the model use fixed-point boundary conditions rather
824 : ! than zero-flux boundary conditions. This means that the value of
825 : ! var_zm stays the same over the course of the timestep at the lower
826 : ! boundary, as well as at the upper boundary.
827 : !
828 : ! In order to discretize the boundary conditions for equations requiring
829 : ! fixed-point boundary conditions, either:
830 : ! a) in the parent subroutine or function (that calls this function),
831 : ! loop over all vertical levels from the second-lowest to the
832 : ! second-highest, ignoring the boundary levels. Then set the values
833 : ! at the boundary levels in the parent subroutine; or
834 : ! b) in the parent subroutine or function, loop over all vertical levels
835 : ! and then overwrite the results at the boundary levels.
836 : !
837 : ! Either way, at the boundary levels, an array with a value of 1 at the
838 : ! main diagonal on the left-hand side and with values of 0 at all other
839 : ! diagonals on the left-hand side will preserve the right-hand side value
840 : ! at that level, thus satisfying the fixed-point boundary conditions.
841 : !
842 : !
843 : ! Conservation Properties:
844 : !
845 : ! When zero-flux boundary conditions are used, this technique of
846 : ! discretizing the eddy diffusion term leads to conservative differencing.
847 : ! When conservative differencing is in place, the column totals for each
848 : ! column in the left-hand side matrix (for the eddy diffusion term) should
849 : ! be equal to 0. This ensures that the total amount of the quantity var_zm
850 : ! over the entire vertical domain is being conserved, meaning that nothing
851 : ! is lost due to diffusional effects.
852 : !
853 : ! To see that this conservation law is satisfied, compute the eddy diffusion
854 : ! of var_zm and integrate vertically. In discretized matrix notation (where
855 : ! "i" stands for the matrix column and "j" stands for the matrix row):
856 : !
857 : ! 0 = Sum_j Sum_i ( 1/invrs_dzm )_i
858 : ! ( invrs_dzm * ((K_zt+nu)*invrs_dzt) )_ij (var_zm)_j.
859 : !
860 : ! The left-hand side matrix, ( invrs_dzm * ((K_zt+nu)*invrs_dzt) )_ij, is
861 : ! partially written below. The sum over i in the above equation removes
862 : ! invrs_dzm everywhere from the matrix below. The sum over j leaves the
863 : ! column totals that are desired.
864 : !
865 : ! Left-hand side matrix contributions from eddy diffusion term; first four
866 : ! vertical levels:
867 : !
868 : ! ---------------------------------------------------------------------->
869 : !k=1 | +invrs_dzm(k) -invrs_dzm(k) 0
870 : ! | *(K_zt(k+1)+nu) *(K_zt(k+1)+nu)
871 : ! | *invrs_dzt(k+1) *invrs_dzt(k+1)
872 : ! |
873 : !k=2 | -invrs_dzm(k) +invrs_dzm(k) -invrs_dzm(k)
874 : ! | *(K_zt(k)+nu) *[ (K_zt(k+1)+nu) *(K_zt(k+1)+nu)
875 : ! | *invrs_dzt(k) *invrs_dzt(k+1) *invrs_dzt(k+1)
876 : ! | +(K_zt(k)+nu)
877 : ! | *invrs_dzt(k) ]
878 : ! |
879 : !k=3 | 0 -invrs_dzm(k) +invrs_dzm(k)
880 : ! | *(K_zt(k)+nu) *[ (K_zt(k+1)+nu)
881 : ! | *invrs_dzt(k) *invrs_dzt(k+1)
882 : ! | +(K_zt(k)+nu)
883 : ! | *invrs_dzt(k) ]
884 : ! |
885 : !k=4 | 0 0 -invrs_dzm(k)
886 : ! | *(K_zt(k)+nu)
887 : ! | *invrs_dzt(k)
888 : ! \ /
889 : !
890 : ! Note: The superdiagonal term from level 3 and both the main diagonal and
891 : ! superdiagonal terms from level 4 are not shown on this diagram.
892 : !
893 : ! Note: The matrix shown is a tridiagonal matrix. For a band diagonal
894 : ! matrix (with 5 diagonals), there would be an extra row between each
895 : ! of the rows shown and an extra column between each of the columns
896 : ! shown. However, for the purposes of the var_zm eddy diffusion
897 : ! term, those extra row and column values are all 0, and the
898 : ! conservation properties of the matrix aren't effected.
899 : !
900 : ! If fixed-point boundary conditions are used, the matrix entries at
901 : ! level 1 (k=1) read: 1 0 0; which means that conservative differencing
902 : ! is not in play. The total amount of var_zm over the entire vertical
903 : ! domain is not being conserved, as amounts of var_zm may be fluxed out
904 : ! through the upper boundary or lower boundary through the effects of
905 : ! diffusion.
906 : !
907 : ! Brian Griffin. April 26, 2008.
908 :
909 : ! References:
910 : ! None
911 : !-----------------------------------------------------------------------
912 :
913 : use grid_class, only: &
914 : grid, & ! Type
915 : ddzt ! Procedure
916 :
917 : use constants_clubb, only: &
918 : zero ! Constant(s)
919 :
920 : use clubb_precision, only: &
921 : core_rknd ! Variable(s)
922 :
923 : use model_flags, only: &
924 : l_upwind_Kh_dp_term
925 :
926 : implicit none
927 :
928 : ! Constant parameters
929 : integer, parameter :: &
930 : kp1_mdiag = 1, & ! Momentum superdiagonal index.
931 : k_mdiag = 2, & ! Momentum main diagonal index.
932 : km1_mdiag = 3 ! Momentum subdiagonal index.
933 :
934 : ! Input Variables
935 : integer, intent(in) :: &
936 : nz, &
937 : ngrdcol
938 :
939 : type (grid), target, intent(in) :: gr
940 :
941 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
942 : K_zm, &
943 : K_zt, & ! Coef. of eddy diffusivity at thermo. levels [m^2/s]
944 : rho_ds_zt, &
945 : invrs_rho_ds_zm
946 :
947 : real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
948 : nu ! Background constant coef. of eddy diffusivity [m^2/s]
949 :
950 : ! Return Variable
951 : real( kind = core_rknd ), dimension(3,ngrdcol,nz), intent(out) :: &
952 : lhs ! LHS coefficient of diffusion term [1/s]
953 :
954 : ! Local Variable
955 : integer :: i, k ! Vertical level index
956 :
957 : real( kind = core_rknd ), dimension(3,ngrdcol,nz) :: &
958 2823552 : lhs_upwind ! LHS coefficient diffusion term due to upwind [1/s]
959 :
960 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
961 2823552 : drhoKdz_zm, &
962 2823552 : rho_K_zt_nu, &
963 2823552 : ddzt_rho_K_zt_nu
964 :
965 : !------------Begin code------------------------------
966 :
967 : !$acc data copyin( gr, gr%invrs_dzt, gr%invrs_dzm, &
968 : !$acc K_zm, K_zt, rho_ds_zt, invrs_rho_ds_zm, nu ) &
969 : !$acc copyout( lhs ) &
970 : !$acc create( lhs_upwind, drhoKdz_zm, rho_K_zt_nu, ddzt_rho_K_zt_nu )
971 :
972 : !$acc parallel loop gang vector collapse(2) default(present)
973 121412736 : do k = 1, nz
974 2005148736 : do i = 1, ngrdcol
975 2003736960 : rho_K_zt_nu(i,k) = rho_ds_zt(i,k) * ( K_zt(i,k) + nu(i) )
976 : end do
977 : end do
978 : !$acc end parallel loop
979 :
980 : ! calculate the dKh_zm/dz
981 1411776 : ddzt_rho_K_zt_nu = ddzt( nz, ngrdcol, gr, rho_K_zt_nu )
982 :
983 : !$acc parallel loop gang vector collapse(2) default(present)
984 121412736 : do k = 1, nz
985 2005148736 : do i = 1, ngrdcol
986 2003736960 : drhoKdz_zm(i,k) = - invrs_rho_ds_zm(i,k) * ddzt_rho_K_zt_nu(i,k)
987 : end do
988 : end do
989 : !$acc end parallel loop
990 :
991 : ! extra terms with upwind scheme
992 : ! k = 1 (bottom level); lowere boundary level
993 : !$acc parallel loop gang vector default(present)
994 23573376 : do i = 1, ngrdcol
995 22161600 : lhs_upwind(kp1_mdiag,i,1) = + min( drhoKdz_zm(i,1) , zero ) * gr%invrs_dzt(i,2)
996 22161600 : lhs_upwind(k_mdiag,i,1) = - min( drhoKdz_zm(i,1) , zero ) * gr%invrs_dzt(i,2)
997 23573376 : lhs_upwind(km1_mdiag,i,1) = zero
998 : end do
999 : !$acc end parallel loop
1000 :
1001 : ! Most of the interior model; normal conditions.
1002 : !$acc parallel loop gang vector collapse(2) default(present)
1003 118589184 : do k = 2, nz-1, 1
1004 1958001984 : do i = 1, ngrdcol
1005 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1006 1839412800 : lhs_upwind(kp1_mdiag,i,k) = + min( drhoKdz_zm(i,k) , zero ) * gr%invrs_dzt(i,k+1)
1007 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1008 0 : lhs_upwind(k_mdiag,i,k) = - min( drhoKdz_zm(i,k) , zero ) * gr%invrs_dzt(i,k+1) &
1009 1839412800 : + max( drhoKdz_zm(i,k) , zero ) * gr%invrs_dzt(i,k)
1010 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1011 1956590208 : lhs_upwind(km1_mdiag,i,k) = - max( drhoKdz_zm(i,k) , zero ) * gr%invrs_dzt(i,k)
1012 : end do
1013 : end do
1014 : !$acc end parallel loop
1015 :
1016 : ! k = nz (top level); upper boundary level.
1017 : ! Only relevant if zero-flux boundary conditions are used.
1018 : !$acc parallel loop gang vector default(present)
1019 23573376 : do i = 1, ngrdcol
1020 22161600 : lhs_upwind(kp1_mdiag,i,nz) = zero
1021 22161600 : lhs_upwind(k_mdiag,i,nz) = + max( drhoKdz_zm(i,nz) , zero ) * gr%invrs_dzt(i,nz)
1022 23573376 : lhs_upwind(km1_mdiag,i,nz) = - max( drhoKdz_zm(i,nz) , zero ) * gr%invrs_dzt(i,nz)
1023 : end do
1024 : !$acc end parallel loop
1025 :
1026 : ! k = 1; lower boundary level at surface.
1027 : ! Only relevant if zero-flux boundary conditions are used.
1028 : ! These k=1 lines currently do not have any effect on model results.
1029 : ! These k=1 level of this "lhs" array is not fed into
1030 : ! the final LHS matrix that will be used to solve for the next timestep.
1031 :
1032 : if ( .not. l_upwind_Kh_dp_term ) then
1033 :
1034 : !$acc parallel loop gang vector default(present)
1035 23573376 : do i = 1, ngrdcol
1036 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1037 44323200 : lhs(kp1_mdiag,i,1) = - gr%invrs_dzm(i,1) * invrs_rho_ds_zm(i,1) &
1038 44323200 : * ( K_zt(i,2) + nu(i) ) * rho_ds_zt(i,2) * gr%invrs_dzt(i,2)
1039 :
1040 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1041 0 : lhs(k_mdiag,i,1) = + gr%invrs_dzm(i,1) * invrs_rho_ds_zm(i,1) &
1042 22161600 : * ( K_zt(i,2) + nu(i) ) * rho_ds_zt(i,2) * gr%invrs_dzt(i,2)
1043 :
1044 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1045 23573376 : lhs(km1_mdiag,i,1) = zero
1046 : end do
1047 : !$acc end parallel loop
1048 :
1049 : else
1050 :
1051 : !$acc parallel loop gang vector default(present)
1052 : do i = 1, ngrdcol
1053 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1054 : lhs(kp1_mdiag,i,1) = - gr%invrs_dzm(i,1) * ( K_zm(i,1) + nu(i) ) * gr%invrs_dzt(i,2) &
1055 : + lhs_upwind(kp1_mdiag,i,1)
1056 :
1057 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1058 : lhs(k_mdiag,i,1) = + gr%invrs_dzm(i,1) * ( K_zm(i,1) + nu(i) ) * gr%invrs_dzt(i,2) &
1059 : + lhs_upwind(k_mdiag,i,1)
1060 :
1061 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1062 : lhs(km1_mdiag,i,1) = zero + lhs_upwind(km1_mdiag,i,1)
1063 : end do
1064 : !$acc end parallel loop
1065 :
1066 : end if
1067 :
1068 : ! Most of the interior model; normal conditions.
1069 : if ( .not. l_upwind_Kh_dp_term ) then
1070 :
1071 : !$acc parallel loop gang vector collapse(2) default(present)
1072 118589184 : do k = 2, nz-1, 1
1073 1958001984 : do i = 1, ngrdcol
1074 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1075 3678825600 : lhs(kp1_mdiag,i,k) &
1076 0 : = - gr%invrs_dzm(i,k) * invrs_rho_ds_zm(i,k) &
1077 5518238400 : * ( K_zt(i,k+1) + nu(i) ) * rho_ds_zt(i,k+1) * gr%invrs_dzt(i,k+1)
1078 :
1079 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1080 : lhs(k_mdiag,i,k) &
1081 0 : = + gr%invrs_dzm(i,k) * invrs_rho_ds_zm(i,k) &
1082 0 : * ( ( K_zt(i,k+1) + nu(i) ) * rho_ds_zt(i,k+1) * gr%invrs_dzt(i,k+1) &
1083 1839412800 : + ( K_zt(i,k) + nu(i) ) * rho_ds_zt(i,k) * gr%invrs_dzt(i,k) )
1084 :
1085 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1086 : lhs(km1_mdiag,i,k) &
1087 0 : = - gr%invrs_dzm(i,k) * invrs_rho_ds_zm(i,k) &
1088 1956590208 : * ( K_zt(i,k) + nu(i) ) * rho_ds_zt(i,k) * gr%invrs_dzt(i,k)
1089 : end do
1090 : end do
1091 : !$acc end parallel loop
1092 :
1093 : else
1094 :
1095 : !$acc parallel loop gang vector collapse(2) default(present)
1096 : do k = 2, nz-1, 1
1097 : do i = 1, ngrdcol
1098 :
1099 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1100 : lhs(kp1_mdiag,i,k) &
1101 : = - gr%invrs_dzm(i,k) * ( K_zm(i,k) + nu(i) ) * gr%invrs_dzt(i,k+1) &
1102 : + lhs_upwind(kp1_mdiag,i,k)
1103 :
1104 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1105 : lhs(k_mdiag,i,k) &
1106 : = + gr%invrs_dzm(i,k) * ( K_zm(i,k) + nu(i) ) * ( gr%invrs_dzt(i,k+1) + gr%invrs_dzt(i,k) ) &
1107 : + lhs_upwind(k_mdiag,i,k)
1108 :
1109 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1110 : lhs(km1_mdiag,i,k) &
1111 : = - gr%invrs_dzm(i,k) * ( K_zm(i,k) + nu(i) ) * gr%invrs_dzt(i,k) &
1112 : + lhs_upwind(km1_mdiag,i,k)
1113 :
1114 : end do
1115 : end do
1116 : !$acc end parallel loop
1117 :
1118 : end if
1119 :
1120 : ! k = nz (top level); upper boundary level.
1121 : ! Only relevant if zero-flux boundary conditions are used.
1122 :
1123 : if ( .not. l_upwind_Kh_dp_term ) then
1124 :
1125 : !$acc parallel loop gang vector default(present)
1126 23573376 : do i = 1, ngrdcol
1127 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1128 22161600 : lhs(kp1_mdiag,i,nz) = zero
1129 :
1130 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1131 : lhs(k_mdiag,i,nz) &
1132 0 : = + gr%invrs_dzm(i,nz) * invrs_rho_ds_zm(i,nz) &
1133 22161600 : * ( K_zt(i,nz) + nu(i) ) * rho_ds_zt(i,nz) * gr%invrs_dzt(i,nz)
1134 :
1135 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1136 : lhs(km1_mdiag,i,nz) &
1137 0 : = - gr%invrs_dzm(i,nz) * invrs_rho_ds_zm(i,nz) &
1138 23573376 : * ( K_zt(i,nz) + nu(i) ) * rho_ds_zt(i,nz) * gr%invrs_dzt(i,nz)
1139 : end do
1140 : !$acc end parallel loop
1141 :
1142 : else
1143 :
1144 : !$acc parallel loop gang vector default(present)
1145 : do i = 1, ngrdcol
1146 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
1147 : lhs(kp1_mdiag,i,nz) = zero &
1148 : + lhs_upwind(kp1_mdiag,i,nz)
1149 :
1150 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
1151 : lhs(k_mdiag,i,nz) &
1152 : = + gr%invrs_dzm(i,nz) * ( K_zm(i,nz) + nu(i) ) * gr%invrs_dzt(i,nz) &
1153 : + lhs_upwind(k_mdiag,i,nz)
1154 :
1155 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
1156 : lhs(km1_mdiag,i,nz) &
1157 : = - gr%invrs_dzm(i,nz) * ( K_zm(i,nz) + nu(i) ) * gr%invrs_dzt(i,nz) &
1158 : + lhs_upwind(km1_mdiag,i,nz)
1159 : end do
1160 : !$acc end parallel loop
1161 :
1162 : end if
1163 :
1164 : !$acc end data
1165 :
1166 1411776 : return
1167 :
1168 : end subroutine diffusion_zm_lhs
1169 :
1170 : !=============================================================================
1171 :
1172 : end module diffusion
|