Line data Source code
1 : !-----------------------------------------------------------------------
2 : ! $Id$
3 : !===============================================================================
4 : module mean_adv
5 :
6 : ! Description:
7 : ! Module mean_adv computes the mean advection terms for all of the
8 : ! time-tendency (prognostic) equations in the CLUBB parameterization. All of
9 : ! the mean advection terms are solved for completely implicitly, and therefore
10 : ! become part of the left-hand side of their respective equations.
11 : !
12 : ! Function term_ma_zt_lhs handles the mean advection terms for the variables
13 : ! located at thermodynamic grid levels. These variables are: rtm, thlm, wp3,
14 : ! all hydrometeor species, and sclrm.
15 : !
16 : ! Function term_ma_zm_lhs handles the mean advection terms for the variables
17 : ! located at momentum grid levels. The variables are: wprtp, wpthlp, wp2,
18 : ! rtp2, thlp2, rtpthlp, up2, vp2, wpsclrp, sclrprtp, sclrpthlp, and sclrp2.
19 :
20 : implicit none
21 :
22 : private ! Default scope
23 :
24 : public :: term_ma_zt_lhs, &
25 : term_ma_zm_lhs
26 :
27 : integer, parameter :: &
28 : ndiags3 = 3
29 :
30 : contains
31 :
32 : !=============================================================================
33 352944 : subroutine term_ma_zt_lhs( nz, ngrdcol, wm_zt, weights_zt2zm, & ! Intent(in)
34 352944 : invrs_dzt, invrs_dzm, & ! Intent(in)
35 : l_upwind_xm_ma, & ! Intent(in)
36 352944 : lhs_ma ) ! Intent(out)
37 :
38 : ! Description:
39 : ! Mean advection 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 a mean advection term:
45 : !
46 : ! - w * d(var_zt)/dz.
47 : !
48 : ! This term is solved for completely implicitly, such that:
49 : !
50 : ! - w * d( var_zt(t+1) )/dz.
51 : !
52 : ! Note: When the term is brought over to the left-hand side, the sign
53 : ! is reversed and the leading "-" in front of the term is changed to
54 : ! a "+".
55 : !
56 : ! The timestep index (t+1) means that the value of var_zt being used is from
57 : ! the next timestep, which is being advanced to in solving the d(var_zt)/dt
58 : ! equation.
59 : !
60 : ! This term is discretized as follows:
61 : !
62 : ! The values of var_zt are found on the thermodynamic levels, as are the
63 : ! values of wm_zt (mean vertical velocity on thermodynamic levels). The
64 : ! variable var_zt is interpolated to the intermediate momentum levels. The
65 : ! derivative of the interpolated values is taken over the central
66 : ! thermodynamic level. The derivative is multiplied by wm_zt at the central
67 : ! thermodynamic level to get the desired result.
68 : !
69 : ! -----var_zt---------------------------------------------- t(k+1)
70 : !
71 : ! =================var_zt(interp)========================== m(k)
72 : !
73 : ! -----var_zt---------------------d(var_zt)/dz-----wm_zt--- t(k)
74 : !
75 : ! =================var_zt(interp)========================== m(k-1)
76 : !
77 : ! -----var_zt---------------------------------------------- t(k-1)
78 : !
79 : ! The vertical indices t(k+1), m(k), t(k), m(k-1), and t(k-1) correspond
80 : ! with altitudes zt(k+1), zm(k), zt(k), zm(k-1), and zt(k-1), respectively.
81 : ! The letter "t" is used for thermodynamic levels and the letter "m" is used
82 : ! for momentum levels.
83 : !
84 : ! invrs_dzt(k) = 1 / ( zm(k) - zm(k-1) )
85 : !
86 : !
87 : ! Special discretization for upper boundary level:
88 : !
89 : ! Method 1: Constant derivative method (or "one-sided" method).
90 : !
91 : ! The values of var_zt are found on the thermodynamic levels, as are the
92 : ! values of wm_zt (mean vertical velocity on the thermodynamic levels). The
93 : ! variable var_zt is interpolated to momentum level gr%nz-1, based on
94 : ! the values of var_zt at thermodynamic levels gr%nz and gr%nz-1.
95 : ! However, the variable var_zt cannot be interpolated to momentum level
96 : ! gr%nz. Rather, a linear extension is used to find the value of var_zt
97 : ! at momentum level gr%nz, based on the values of var_zt at thermodynamic
98 : ! levels gr%nz and gr%nz-1. The derivative of the extended and
99 : ! interpolated values, d(var_zt)/dz, is taken over the central thermodynamic
100 : ! level. Of course, this derivative will be the same as the derivative of
101 : ! var_zt between thermodynamic levels gr%nz and gr%nz-1. The derivative
102 : ! is multiplied by wm_zt at the central thermodynamic level to get the
103 : ! desired result.
104 : !
105 : ! For the following diagram, k = gr%nz, which is the uppermost level of
106 : ! the model:
107 : !
108 : ! =================var_zt(extend)========================== m(k) Boundary
109 : !
110 : ! -----var_zt---------------------d(var_zt)/dz-----wm_zt--- t(k)
111 : !
112 : ! =================var_zt(interp)========================== m(k-1)
113 : !
114 : ! -----var_zt---------------------------------------------- t(k-1)
115 : !
116 : !
117 : ! Method 2: Zero derivative method:
118 : ! the derivative d(var_zt)/dz over the model top is set to 0.
119 : !
120 : ! This method corresponds with the "zero-flux" boundary condition option
121 : ! for eddy diffusion, where d(var_zt)/dz is set to 0 across the upper
122 : ! boundary.
123 : !
124 : ! In order to discretize the upper boundary condition, consider a new level
125 : ! outside the model (thermodynamic level gr%nz+1) just above the upper
126 : ! boundary level (thermodynamic level gr%nz). The value of var_zt at the
127 : ! level just outside the model is defined to be the same as the value of
128 : ! var_zt at thermodynamic level gr%nz. Therefore, the value of
129 : ! d(var_zt)/dz between the level just outside the model and the uppermost
130 : ! thermodynamic level is 0, staying consistent with the zero-flux boundary
131 : ! condition option for the eddy diffusion portion of the code. Therefore,
132 : ! the value of var_zt at momentum level gr%nz, which is the upper boundary
133 : ! of the model, would be the same as the value of var_zt at the uppermost
134 : ! thermodynamic level.
135 : !
136 : ! The values of var_zt are found on the thermodynamic levels, as are the
137 : ! values of wm_zt (mean vertical velocity on the thermodynamic levels). The
138 : ! variable var_zt is interpolated to momentum level gr%nz-1, based on
139 : ! the values of var_zt at thermodynamic levels gr%nz and gr%nz-1. The
140 : ! value of var_zt at momentum level gr%nz is set equal to the value of
141 : ! var_zt at thermodynamic level gr%nz, as described above. The derivative
142 : ! of the set and interpolated values, d(var_zt)/dz, is taken over the
143 : ! central thermodynamic level. The derivative is multiplied by wm_zt at the
144 : ! central thermodynamic level to get the desired result.
145 : !
146 : ! For the following diagram, k = gr%nz, which is the uppermost level of
147 : ! the model:
148 : !
149 : ! --[var_zt(k+1) = var_zt(k)]----(level outside model)----- t(k+1)
150 : !
151 : ! ==[var_zt(top) = var_zt(k)]===[d(var_zt)/dz|_(top) = 0]== m(k) Boundary
152 : !
153 : ! -----var_zt(k)------------------d(var_zt)/dz-----wm_zt--- t(k)
154 : !
155 : ! =================var_zt(interp)========================== m(k-1)
156 : !
157 : ! -----var_zt(k-1)----------------------------------------- t(k-1)
158 : !
159 : ! where (top) stands for the grid index of momentum level k = gr%nz, which
160 : ! is the upper boundary of the model.
161 : !
162 : ! This method of boundary discretization is also similar to the method
163 : ! currently employed at the lower boundary for most thermodynamic-level
164 : ! variables. Since thermodynamic level k = 1 is below the model bottom,
165 : ! mean advection is not applied. Thus, thermodynamic level k = 2 becomes
166 : ! the lower boundary level. Now, the mean advection term at thermodynamic
167 : ! level 2 takes into account var_zt from levels 1, 2, and 3. However, in
168 : ! most cases, the value of var_zt(1) is set equal to var_zt(2) after the
169 : ! matrix of equations has been solved. Therefore, the derivative,
170 : ! d(var_zt)/dz, over the model bottom (momentum level k = 1) becomes 0.
171 : ! Thus, the method of setting d(var_zt)/dz to 0 over the model top keeps
172 : ! the way the upper and lower boundaries are handled consistent with each
173 : ! other.
174 :
175 : ! References:
176 : ! None
177 : !-----------------------------------------------------------------------
178 :
179 : use grid_class, only: &
180 : grid ! Type
181 :
182 : use constants_clubb, only: &
183 : one, & ! Constant(s)
184 : zero
185 :
186 : use clubb_precision, only: &
187 : core_rknd ! Variable(s)
188 :
189 : implicit none
190 :
191 : ! Constant parameters
192 : integer, parameter :: &
193 : kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
194 : k_tdiag = 2, & ! Thermodynamic main diagonal index.
195 : km1_tdiag = 3 ! Thermodynamic subdiagonal index.
196 :
197 : integer, parameter :: &
198 : t_above = 1, & ! Index for upper thermodynamic level grid weight.
199 : t_below = 2 ! Index for lower thermodynamic level grid weight.
200 :
201 : ! -------------------------- Input Variables --------------------------
202 : integer, intent(in) :: &
203 : nz, &
204 : ngrdcol
205 :
206 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
207 : wm_zt, & ! wm_zt [m/s]
208 : invrs_dzt, & ! Inverse of grid spacing [1/m]
209 : invrs_dzm ! Inverse of grid spacing [1/m]
210 :
211 : real( kind = core_rknd ), dimension(ngrdcol,nz,t_above:t_below), intent(in) :: &
212 : weights_zt2zm
213 :
214 : logical, intent(in) :: &
215 : l_upwind_xm_ma ! This flag determines whether we want to use an upwind
216 : ! differencing approximation rather than a centered
217 : ! differencing for turbulent or mean advection terms.
218 : ! It affects rtm, thlm, sclrm, um and vm.
219 :
220 : ! -------------------------- Return Variable --------------------------
221 : real( kind = core_rknd ), dimension(ndiags3,ngrdcol,nz), intent(out) :: &
222 : lhs_ma ! Mean advection contributions to lhs [1/s]
223 :
224 : ! -------------------------- Local Variables --------------------------
225 :
226 : integer :: i, k, b ! Vertical level index
227 :
228 : !-------------------------- Begin Code --------------------------
229 :
230 : !$acc data copyin( wm_zt, invrs_dzt, invrs_dzm, weights_zt2zm ) &
231 : !$acc copyout( lhs_ma )
232 :
233 : ! Set lower boundary array to 0
234 : !$acc parallel loop gang vector collapse(2) default(present)
235 5893344 : do i = 1, ngrdcol
236 22514544 : do b = 1, ndiags3
237 22161600 : lhs_ma(b,i,1) = 0.0_core_rknd
238 : end do
239 : end do
240 : !$acc end parallel loop
241 :
242 :
243 352944 : if ( .not. l_upwind_xm_ma ) then ! Use centered differencing
244 :
245 : ! Most of the interior model; normal conditions.
246 : !$acc parallel loop gang vector collapse(2) default(present)
247 0 : do k = 2, nz, 1
248 0 : do i = 1, ngrdcol
249 :
250 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
251 0 : lhs_ma(kp1_tdiag,i,k) = + wm_zt(i,k) * invrs_dzt(i,k) * weights_zt2zm(i,k,t_above)
252 :
253 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
254 : lhs_ma(k_tdiag,i,k) = + wm_zt(i,k) * invrs_dzt(i,k) * ( weights_zt2zm(i,k,t_below) &
255 0 : - weights_zt2zm(i,k-1,t_above) )
256 :
257 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
258 0 : lhs_ma(km1_tdiag,i,k) = - wm_zt(i,k) * invrs_dzt(i,k) * weights_zt2zm(i,k-1,t_below)
259 : end do
260 : end do ! k = 2, nz, 1
261 : !$acc end parallel loop
262 :
263 : ! Upper Boundary
264 :
265 : ! Special discretization for zero derivative method, where the
266 : ! derivative d(var_zt)/dz over the model top is set to 0, in order
267 : ! to stay consistent with the zero-flux boundary condition option
268 : ! in the eddy diffusion code.
269 : !$acc parallel loop gang vector default(present)
270 0 : do i = 1, ngrdcol
271 :
272 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
273 0 : lhs_ma(kp1_tdiag,i,nz) = zero
274 :
275 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
276 : lhs_ma(k_tdiag,i,nz) = + wm_zt(i,nz) * invrs_dzt(i,nz) &
277 0 : * ( one - weights_zt2zm(i,nz-1,t_above) )
278 :
279 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
280 : lhs_ma(km1_tdiag,i,nz) = - wm_zt(i,nz) * invrs_dzt(i,nz) &
281 0 : * weights_zt2zm(i,nz-1,t_below)
282 : end do
283 : !$acc end parallel loop
284 :
285 : else ! l_upwind_xm_ma == .true.; use "upwind" differencing
286 :
287 : ! Most of the interior model; normal conditions.
288 : !$acc parallel loop gang vector collapse(2) default(present)
289 30000240 : do k = 2, nz, 1
290 495393840 : do i = 1, ngrdcol
291 495040896 : if ( wm_zt(i,k) >= zero ) then ! Mean wind is in upward direction
292 :
293 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
294 212765418 : lhs_ma(kp1_tdiag,i,k) = zero
295 :
296 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
297 212765418 : lhs_ma(k_tdiag,i,k) = + wm_zt(i,k) * invrs_dzm(i,k-1)
298 :
299 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
300 212765418 : lhs_ma(km1_tdiag,i,k) = - wm_zt(i,k) * invrs_dzm(i,k-1)
301 :
302 : else ! wm_zt < 0; Mean wind is in downward direction
303 :
304 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
305 252628182 : lhs_ma(kp1_tdiag,i,k) = + wm_zt(i,k) * invrs_dzm(i,k)
306 :
307 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
308 252628182 : lhs_ma(k_tdiag,i,k) = - wm_zt(i,k) * invrs_dzm(i,k)
309 :
310 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
311 252628182 : lhs_ma(km1_tdiag,i,k) = zero
312 :
313 : endif ! wm_zt > 0
314 :
315 : end do
316 : end do ! k = 2, nz, 1
317 : !$acc end parallel loop
318 :
319 : ! Upper Boundary
320 : !$acc parallel loop gang vector default(present)
321 5893344 : do i = 1, ngrdcol
322 5893344 : if ( wm_zt(i,nz) >= zero ) then ! Mean wind is in upward direction
323 :
324 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
325 2608740 : lhs_ma(kp1_tdiag,i,nz) = zero
326 :
327 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
328 2608740 : lhs_ma(k_tdiag,i,nz) = + wm_zt(i,nz) * invrs_dzm(i,nz-1)
329 :
330 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
331 2608740 : lhs_ma(km1_tdiag,i,nz) = - wm_zt(i,nz) * invrs_dzm(i,nz-1)
332 :
333 : else ! wm_zt < 0; Mean wind is in downward direction
334 :
335 : ! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
336 2931660 : lhs_ma(kp1_tdiag,i,nz) = zero
337 :
338 : ! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
339 2931660 : lhs_ma(k_tdiag,i,nz) = zero
340 :
341 : ! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
342 2931660 : lhs_ma(km1_tdiag,i,nz) = zero
343 :
344 : end if ! wm_zt > 0
345 : end do
346 : !$acc end parallel loop
347 :
348 : endif ! l_upwind_xm_ma
349 :
350 : !$acc end data
351 :
352 352944 : return
353 :
354 : end subroutine term_ma_zt_lhs
355 :
356 : !=============================================================================
357 1058832 : subroutine term_ma_zm_lhs( nz, ngrdcol, wm_zm, &
358 1058832 : invrs_dzm, weights_zm2zt, &
359 1058832 : lhs_ma )
360 :
361 : ! Description:
362 : ! Mean advection of var_zm: implicit portion of the code.
363 : !
364 : ! The variable "var_zm" stands for a variable that is located at momentum
365 : ! grid levels.
366 : !
367 : ! The d(var_zm)/dt equation contains a mean advection term:
368 : !
369 : ! - w * d(var_zm)/dz.
370 : !
371 : ! This term is solved for completely implicitly, such that:
372 : !
373 : ! - w * d( var_zm(t+1) )/dz.
374 : !
375 : ! Note: When the term is brought over to the left-hand side, the sign
376 : ! is reversed and the leading "-" in front of the term is changed to
377 : ! a "+".
378 : !
379 : ! The timestep index (t+1) means that the value of var_zm being used is from
380 : ! the next timestep, which is being advanced to in solving the d(var_zm)/dt
381 : ! equation.
382 : !
383 : ! This term is discretized as follows:
384 : !
385 : ! The values of var_zm are found on the momentum levels, as are the values
386 : ! of wm_zm (mean vertical velocity on momentum levels). The variable var_zm
387 : ! is interpolated to the intermediate thermodynamic levels. The derivative
388 : ! of the interpolated values is taken over the central momentum level. The
389 : ! derivative is multiplied by wm_zm at the central momentum level to get the
390 : ! desired result.
391 : !
392 : ! =====var_zm============================================== m(k+1)
393 : !
394 : ! -----------------var_zm(interp)-------------------------- t(k+1)
395 : !
396 : ! =====var_zm=====================d(var_zm)/dz=====wm_zm=== m(k)
397 : !
398 : ! -----------------var_zm(interp)-------------------------- t(k)
399 : !
400 : ! =====var_zm============================================== m(k-1)
401 : !
402 : ! The vertical indices m(k+1), t(k+1), m(k), t(k), and m(k-1) correspond
403 : ! with altitudes zm(k+1), zt(k+1), zm(k), zt(k), and zm(k-1), respectively.
404 : ! The letter "t" is used for thermodynamic levels and the letter "m" is used
405 : ! for momentum levels.
406 : !
407 : ! invrs_dzm(k) = 1 / ( zt(k+1) - zt(k) )
408 :
409 : ! References:
410 : !-----------------------------------------------------------------------
411 :
412 : use grid_class, only: &
413 : grid ! Type
414 :
415 : use constants_clubb, only: &
416 : zero ! Constant(s)
417 :
418 : use clubb_precision, only: &
419 : core_rknd ! Variable(s)
420 :
421 : implicit none
422 :
423 : ! Constant parameters
424 : integer, parameter :: &
425 : kp1_mdiag = 1, & ! Momentum superdiagonal index.
426 : k_mdiag = 2, & ! Momentum main diagonal index.
427 : km1_mdiag = 3 ! Momentum subdiagonal index.
428 :
429 : integer, parameter :: &
430 : m_above = 1, & ! Index for upper momentum level grid weight.
431 : m_below = 2 ! Index for lower momentum level grid weight.
432 :
433 : ! -------------------------- Input Variables --------------------------
434 : integer, intent(in) :: &
435 : nz, &
436 : ngrdcol
437 :
438 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
439 : wm_zm, & ! wm_zm [m/s]
440 : invrs_dzm ! Inverse of grid spacing [1/m]
441 :
442 : real( kind = core_rknd ), dimension(ngrdcol,nz,m_above:m_below), intent(in) :: &
443 : weights_zm2zt
444 :
445 : ! -------------------------- Return Variable --------------------------
446 : real( kind = core_rknd ), dimension(ndiags3,ngrdcol,nz), intent(out) :: &
447 : lhs_ma ! Mean advection contributions to lhs [1/s]
448 :
449 : ! -------------------------- Local Variables
450 : integer :: i, k, b ! Vertical level index
451 :
452 : ! -------------------------- Begin Code --------------------------
453 :
454 : !$acc data copyin( wm_zm, invrs_dzm, weights_zm2zt ) &
455 : !$acc copyout( lhs_ma )
456 :
457 : ! Set lower boundary array to 0
458 : !$acc parallel loop gang vector collapse(2) default(present)
459 17680032 : do i = 1, ngrdcol
460 67543632 : do b = 1, ndiags3
461 66484800 : lhs_ma(b,i,1) = zero
462 : end do
463 : end do
464 : !$acc end parallel loop
465 :
466 : ! Most of the interior model; normal conditions.
467 : !$acc parallel loop gang vector collapse(2) default(present)
468 88941888 : do k = 2, nz-1, 1
469 1468501488 : do i = 1, ngrdcol
470 :
471 : ! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
472 1379559600 : lhs_ma(kp1_mdiag,i,k) = + wm_zm(i,k) * invrs_dzm(i,k) * weights_zm2zt(i,k+1,m_above)
473 :
474 : ! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
475 : lhs_ma(k_mdiag,i,k) = + wm_zm(i,k) * invrs_dzm(i,k) * ( weights_zm2zt(i,k+1,m_below) &
476 1379559600 : - weights_zm2zt(i,k,m_above) )
477 :
478 : ! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
479 1467442656 : lhs_ma(km1_mdiag,i,k) = - wm_zm(i,k) * invrs_dzm(i,k) * weights_zm2zt(i,k,m_below)
480 :
481 : end do
482 : end do ! k = 2, nz-1, 1
483 : !$acc end parallel loop
484 :
485 : ! Set upper boundary array to 0
486 : !$acc parallel loop gang vector collapse(2) default(present)
487 17680032 : do i = 1, ngrdcol
488 67543632 : do b = 1, ndiags3
489 66484800 : lhs_ma(b,i,nz) = zero
490 : end do
491 : end do
492 : !$acc end parallel loop
493 :
494 : !$acc end data
495 :
496 1058832 : return
497 :
498 : end subroutine term_ma_zm_lhs
499 :
500 : !=============================================================================
501 :
502 : end module mean_adv
|