forked from ESCOMP/CLUBB_CESM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
mean_adv.F90
764 lines (565 loc) · 26.5 KB
/
mean_adv.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
!-----------------------------------------------------------------------
! $Id$
!===============================================================================
module mean_adv
! Description:
! Module mean_adv computes the mean advection terms for all of the
! time-tendency (prognostic) equations in the CLUBB parameterization. All of
! the mean advection terms are solved for completely implicitly, and therefore
! become part of the left-hand side of their respective equations.
!
! Function term_ma_zt_lhs handles the mean advection terms for the variables
! located at thermodynamic grid levels. These variables are: rtm, thlm, wp3,
! all hydrometeor species, and sclrm.
!
! Function term_ma_zm_lhs handles the mean advection terms for the variables
! located at momentum grid levels. The variables are: wprtp, wpthlp, wp2,
! rtp2, thlp2, rtpthlp, up2, vp2, wpsclrp, sclrprtp, sclrpthlp, and sclrp2.
implicit none
private ! Default scope
public :: term_ma_zt_lhs, &
term_ma_zt_lhs_all, &
term_ma_zm_lhs, &
term_ma_zm_lhs_all
contains
!=============================================================================
pure function term_ma_zt_lhs( wm_zt, invrs_dzt, level, &
invrs_dzm_k, invrs_dzm_km1 ) &
result( lhs )
! Description:
! Mean advection of var_zt: implicit portion of the code.
!
! The variable "var_zt" stands for a variable that is located at
! thermodynamic grid levels.
!
! The d(var_zt)/dt equation contains a mean advection term:
!
! - w * d(var_zt)/dz.
!
! This term is solved for completely implicitly, such that:
!
! - w * d( var_zt(t+1) )/dz.
!
! Note: When the term is brought over to the left-hand side, the sign
! is reversed and the leading "-" in front of the term is changed to
! a "+".
!
! The timestep index (t+1) means that the value of var_zt being used is from
! the next timestep, which is being advanced to in solving the d(var_zt)/dt
! equation.
!
! This term is discretized as follows:
!
! The values of var_zt are found on the thermodynamic levels, as are the
! values of wm_zt (mean vertical velocity on thermodynamic levels). The
! variable var_zt is interpolated to the intermediate momentum levels. The
! derivative of the interpolated values is taken over the central
! thermodynamic level. The derivative is multiplied by wm_zt at the central
! thermodynamic level to get the desired result.
!
! -----var_zt(kp1)----------------------------------------- t(k+1)
!
! =================var_zt(interp)========================== m(k)
!
! -----var_zt(k)------------------d(var_zt)/dz-----wm_zt--- t(k)
!
! =================var_zt(interp)========================== m(k-1)
!
! -----var_zt(km1)----------------------------------------- t(k-1)
!
! The vertical indices t(k+1), m(k), t(k), m(k-1), and t(k-1) correspond
! with altitudes zt(k+1), zm(k), zt(k), zm(k-1), and zt(k-1), respectively.
! The letter "t" is used for thermodynamic levels and the letter "m" is used
! for momentum levels.
!
! invrs_dzt(k) = 1 / ( zm(k) - zm(k-1) )
!
!
! Special discretization for upper boundary level:
!
! Method 1: Constant derivative method (or "one-sided" method).
!
! The values of var_zt are found on the thermodynamic levels, as are the
! values of wm_zt (mean vertical velocity on the thermodynamic levels). The
! variable var_zt is interpolated to momentum level gr%nz-1, based on
! the values of var_zt at thermodynamic levels gr%nz and gr%nz-1.
! However, the variable var_zt cannot be interpolated to momentum level
! gr%nz. Rather, a linear extension is used to find the value of var_zt
! at momentum level gr%nz, based on the values of var_zt at thermodynamic
! levels gr%nz and gr%nz-1. The derivative of the extended and
! interpolated values, d(var_zt)/dz, is taken over the central thermodynamic
! level. Of course, this derivative will be the same as the derivative of
! var_zt between thermodynamic levels gr%nz and gr%nz-1. The derivative
! is multiplied by wm_zt at the central thermodynamic level to get the
! desired result.
!
! For the following diagram, k = gr%nz, which is the uppermost level of
! the model:
!
! =================var_zt(extend)========================== m(k) Boundary
!
! -----var_zt(k)------------------d(var_zt)/dz-----wm_zt--- t(k)
!
! =================var_zt(interp)========================== m(k-1)
!
! -----var_zt(km1)----------------------------------------- t(k-1)
!
!
! Method 2: Zero derivative method:
! the derivative d(var_zt)/dz over the model top is set to 0.
!
! This method corresponds with the "zero-flux" boundary condition option
! for eddy diffusion, where d(var_zt)/dz is set to 0 across the upper
! boundary.
!
! In order to discretize the upper boundary condition, consider a new level
! outside the model (thermodynamic level gr%nz+1) just above the upper
! boundary level (thermodynamic level gr%nz). The value of var_zt at the
! level just outside the model is defined to be the same as the value of
! var_zt at thermodynamic level gr%nz. Therefore, the value of
! d(var_zt)/dz between the level just outside the model and the uppermost
! thermodynamic level is 0, staying consistent with the zero-flux boundary
! condition option for the eddy diffusion portion of the code. Therefore,
! the value of var_zt at momentum level gr%nz, which is the upper boundary
! of the model, would be the same as the value of var_zt at the uppermost
! thermodynamic level.
!
! The values of var_zt are found on the thermodynamic levels, as are the
! values of wm_zt (mean vertical velocity on the thermodynamic levels). The
! variable var_zt is interpolated to momentum level gr%nz-1, based on
! the values of var_zt at thermodynamic levels gr%nz and gr%nz-1. The
! value of var_zt at momentum level gr%nz is set equal to the value of
! var_zt at thermodynamic level gr%nz, as described above. The derivative
! of the set and interpolated values, d(var_zt)/dz, is taken over the
! central thermodynamic level. The derivative is multiplied by wm_zt at the
! central thermodynamic level to get the desired result.
!
! For the following diagram, k = gr%nz, which is the uppermost level of
! the model:
!
! --[var_zt(kp1) = var_zt(k)]----(level outside model)----- t(k+1)
!
! ==[var_zt(top) = var_zt(k)]===[d(var_zt)/dz|_(top) = 0]== m(k) Boundary
!
! -----var_zt(k)------------------d(var_zt)/dz-----wm_zt--- t(k)
!
! =================var_zt(interp)========================== m(k-1)
!
! -----var_zt(km1)----------------------------------------- t(k-1)
!
! where (top) stands for the grid index of momentum level k = gr%nz, which
! is the upper boundary of the model.
!
! This method of boundary discretization is also similar to the method
! currently employed at the lower boundary for most thermodynamic-level
! variables. Since thermodynamic level k = 1 is below the model bottom,
! mean advection is not applied. Thus, thermodynamic level k = 2 becomes
! the lower boundary level. Now, the mean advection term at thermodynamic
! level 2 takes into account var_zt from levels 1, 2, and 3. However, in
! most cases, the value of var_zt(1) is set equal to var_zt(2) after the
! matrix of equations has been solved. Therefore, the derivative,
! d(var_zt)/dz, over the model bottom (momentum level k = 1) becomes 0.
! Thus, the method of setting d(var_zt)/dz to 0 over the model top keeps
! the way the upper and lower boundaries are handled consistent with each
! other.
! References:
! None
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable(s)
use constants_clubb, only: &
one, & ! Constant(s)
zero
use model_flags, only: &
l_upwind_xm_ma ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Constant parameters
integer, parameter :: &
kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
k_tdiag = 2, & ! Thermodynamic main diagonal index.
km1_tdiag = 3 ! Thermodynamic subdiagonal index.
integer, parameter :: &
t_above = 1, & ! Index for upper thermodynamic level grid weight.
t_below = 2 ! Index for lower thermodynamic level grid weight.
! Input Variables
real( kind = core_rknd ), intent(in) :: &
wm_zt, & ! wm_zt(k) [m/s]
invrs_dzt, & ! Inverse of grid spacing (k) [1/m]
invrs_dzm_k, & ! Inverse of grid spacing (k) [1/m]
invrs_dzm_km1 ! Inverse of grid spacing (k-1) [1/m]
integer, intent(in) :: &
level ! Central thermodynamic level (on which calculation occurs).
! Return Variable
real( kind = core_rknd ), dimension(3) :: lhs
! Local Variables
logical, parameter :: &
l_ub_const_deriv = .true. ! Flag to use the "one-sided" upper boundary.
integer :: &
mk, & ! Momentum level directly above central thermodynamic level.
mkm1 ! Momentum level directly below central thermodynamic level.
! Momentum level (k) is between thermodynamic level (k+1)
! and thermodynamic level (k).
mk = level
! Momentum level (k-1) is between thermodynamic level (k)
! and thermodynamic level (k-1).
mkm1 = level - 1
if ( level == 1 ) then
! k = 1 (bottom level); lower boundary level.
! Thermodynamic level k = 1 is below the model bottom, so all effects
! are shut off.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= zero
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= zero
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= zero
elseif ( level > 1 .and. level < gr%nz ) then
! Most of the interior model; normal conditions.
if( .not. l_upwind_xm_ma ) then ! Use centered differencing
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= + wm_zt * invrs_dzt * gr%weights_zt2zm(t_above,mk)
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= + wm_zt * invrs_dzt * ( gr%weights_zt2zm(t_below,mk) &
- gr%weights_zt2zm(t_above,mkm1) )
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= - wm_zt * invrs_dzt * gr%weights_zt2zm(t_below,mkm1)
else ! l_upwind_xm_ma == .true.; use "upwind" differencing
if ( wm_zt >= zero ) then ! Mean wind is in upward direction
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= zero
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= + wm_zt * invrs_dzm_km1
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= - wm_zt * invrs_dzm_km1
else ! wm_zt < 0; Mean wind is in downward direction
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= + wm_zt * invrs_dzm_k
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= - wm_zt * invrs_dzm_k
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= zero
endif ! wm_zt > 0
endif ! l_upwind_xm_ma
elseif ( level == gr%nz ) then
! k = gr%nz (top level); upper boundary level.
if( .not. l_upwind_xm_ma ) then ! Use "centered" differencing
if ( l_ub_const_deriv ) then
! Special discretization for constant derivative method (or
! "one-sided" derivative method).
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= zero
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= + wm_zt * invrs_dzt * ( gr%weights_zt2zm(t_above,mk) &
- gr%weights_zt2zm(t_above,mkm1) )
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= + wm_zt * invrs_dzt * ( gr%weights_zt2zm(t_below,mk) &
- gr%weights_zt2zm(t_below,mkm1) )
else
! Special discretization for zero derivative method, where the
! derivative d(var_zt)/dz over the model top is set to 0, in order
! to stay consistent with the zero-flux boundary condition option
! in the eddy diffusion code.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= zero
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= + wm_zt * invrs_dzt * ( one - gr%weights_zt2zm(t_above,mkm1) )
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= - wm_zt * invrs_dzt * gr%weights_zt2zm(t_below,mkm1)
endif ! l_ub_const_deriv
else ! l_upwind_xm_ma == .true.; use "upwind" differencing
if ( wm_zt >= zero ) then ! Mean wind is in upward direction
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= zero
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= + wm_zt * invrs_dzm_km1
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= - wm_zt * invrs_dzm_km1
else ! wm_zt < 0; Mean wind is in downward direction
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) &
= zero
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) &
= zero
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) &
= zero
endif ! wm_zt > 0
endif ! l_upwind_xm_ma
endif ! level = gr%nz
return
end function term_ma_zt_lhs
!=============================================================================================
pure subroutine term_ma_zt_lhs_all( wm_zt, invrs_dzt, invrs_dzm, & ! Intent(in)
lhs_ma ) ! Intent(out)
! Description:
! This subroutine is an optimized version of term_ma_zt_lhs. term_ma_zt_lhs
! returns a single 3 dimensional array for any specified grid level. This subroutine returns
! an array of 3 dimensional arrays, one for every grid level not including boundary values.
!
! Notes:
! This subroutine exists for performance concerns. It returns all lhs arrays at once
! so that it can be properly vectorized, see clubb:ticket:834 for detail.
!
!---------------------------------------------------------------------------------------------
use grid_class, only: &
gr ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
use model_flags, only: &
l_upwind_xm_ma ! Variable(s)
implicit none
!------------------- Input Variables -------------------
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
wm_zt, & ! wm_zt(k) [m/s]
invrs_dzt, & ! Inverse of grid spacing (k) [1/m]
invrs_dzm
!------------------- Output Variables -------------------
real( kind = core_rknd ), dimension(3,gr%nz), intent(out) :: &
lhs_ma
!---------------- Local Variables -------------------
integer :: &
k ! Loop variable for current grid level
logical, parameter :: &
l_ub_const_deriv = .true. ! Flag to use the "one-sided" upper boundary.
!---------------- Begin Code -------------------
! Set lower boundary array to 0
lhs_ma(:,1) = 0.0_core_rknd
if( .not. l_upwind_xm_ma ) then ! Use "centered" differencing
! Most of the interior model; normal conditions.
do k = 2, gr%nz-1
lhs_ma(1,k) = + wm_zt(k) * invrs_dzt(k) * gr%weights_zt2zm(1,k)
lhs_ma(2,k) = + wm_zt(k) * invrs_dzt(k) * ( gr%weights_zt2zm(2,k) &
- gr%weights_zt2zm(1,k-1) )
lhs_ma(3,k) = - wm_zt(k) * invrs_dzt(k) * gr%weights_zt2zm(2,k-1)
end do
! Upper Boundary
if ( l_ub_const_deriv ) then
lhs_ma(1,gr%nz) = 0.0_core_rknd
lhs_ma(2,gr%nz) = + wm_zt(k) * invrs_dzt(k) * ( gr%weights_zt2zm(1,k) &
- gr%weights_zt2zm(1,k-1) )
lhs_ma(3,gr%nz) = + wm_zt(k) * invrs_dzt(k) * ( gr%weights_zt2zm(2,k) &
- gr%weights_zt2zm(2,k-1) )
else
lhs_ma(1,gr%nz) = 0.0_core_rknd
lhs_ma(2,gr%nz) = + wm_zt(k) * invrs_dzt(k) &
* ( 1.0_core_rknd - gr%weights_zt2zm(1,k-1) )
lhs_ma(3,gr%nz) = - wm_zt(k) * invrs_dzt(k) * gr%weights_zt2zm(2,k-1)
endif ! l_ub_const_deriv
else ! l_upwind_xm_ma == .true.; use "upwind" differencing
do k = 2, gr%nz-1
if ( wm_zt(k) >= 0.0_core_rknd ) then ! Mean wind is in upward direction
lhs_ma(1,k) = 0.0_core_rknd
lhs_ma(2,k) = + wm_zt(k) * invrs_dzm(k-1)
lhs_ma(3,k) = - wm_zt(k) * invrs_dzm(k-1)
else ! wm_zt < 0; Mean wind is in downward direction
lhs_ma(1,k) = + wm_zt(k) * invrs_dzm(k)
lhs_ma(2,k) = - wm_zt(k) * invrs_dzm(k)
lhs_ma(3,k) = 0.0_core_rknd
endif ! wm_zt > 0
end do
! Upper Boundary
if ( wm_zt(k) >= 0.0_core_rknd ) then ! Mean wind is in upward direction
lhs_ma(1,gr%nz) = 0.0_core_rknd
lhs_ma(2,gr%nz) = + wm_zt(k) * invrs_dzm(k-1)
lhs_ma(3,gr%nz) = - wm_zt(k) * invrs_dzm(k-1)
else ! wm_zt < 0; Mean wind is in downward direction
lhs_ma(1,gr%nz) = 0.0_core_rknd
lhs_ma(2,gr%nz) = 0.0_core_rknd
lhs_ma(3,gr%nz) = 0.0_core_rknd
endif ! wm_zt > 0
endif ! l_upwind_xm_ma
return
end subroutine term_ma_zt_lhs_all
!=============================================================================
pure function term_ma_zm_lhs( wm_zm, invrs_dzm, level ) &
result( lhs )
! Description:
! Mean advection of var_zm: implicit portion of the code.
!
! The variable "var_zm" stands for a variable that is located at momentum
! grid levels.
!
! The d(var_zm)/dt equation contains a mean advection term:
!
! - w * d(var_zm)/dz.
!
! This term is solved for completely implicitly, such that:
!
! - w * d( var_zm(t+1) )/dz.
!
! Note: When the term is brought over to the left-hand side, the sign
! is reversed and the leading "-" in front of the term is changed to
! a "+".
!
! The timestep index (t+1) means that the value of var_zm being used is from
! the next timestep, which is being advanced to in solving the d(var_zm)/dt
! equation.
!
! This term is discretized as follows:
!
! The values of var_zm are found on the momentum levels, as are the values
! of wm_zm (mean vertical velocity on momentum levels). The variable var_zm
! is interpolated to the intermediate thermodynamic levels. The derivative
! of the interpolated values is taken over the central momentum level. The
! derivative is multiplied by wm_zm at the central momentum level to get the
! desired result.
!
! =====var_zm(kp1)========================================= m(k+1)
!
! -----------------var_zm(interp)-------------------------- t(k+1)
!
! =====var_zm(k)==================d(var_zm)/dz=====wm_zm=== m(k)
!
! -----------------var_zm(interp)-------------------------- t(k)
!
! =====var_zm(km1)========================================= m(k-1)
!
! The vertical indices m(k+1), t(k+1), m(k), t(k), and m(k-1) correspond
! with altitudes zm(k+1), zt(k+1), zm(k), zt(k), and zm(k-1), respectively.
! The letter "t" is used for thermodynamic levels and the letter "m" is used
! for momentum levels.
!
! invrs_dzm(k) = 1 / ( zt(k+1) - zt(k) )
! References:
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable(s)
use constants_clubb, only: &
zero ! Constant(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Constant parameters
integer, parameter :: &
kp1_mdiag = 1, & ! Momentum superdiagonal index.
k_mdiag = 2, & ! Momentum main diagonal index.
km1_mdiag = 3 ! Momentum subdiagonal index.
integer, parameter :: &
m_above = 1, & ! Index for upper momentum level grid weight.
m_below = 2 ! Index for lower momentum level grid weight.
! Input Variables
real( kind = core_rknd ), intent(in) :: &
wm_zm, & ! wm_zm(k) [m/s]
invrs_dzm ! Inverse of grid spacing (k) [1/m]
integer, intent(in) :: &
level ! Central momentum level (on which calculation occurs).
! Return Variable
real( kind = core_rknd ), dimension(3) :: lhs
! Local Variables
integer :: &
tkp1, & ! Thermodynamic level directly above central momentum level.
tk ! Thermodynamic level directly below central momentum level.
! Thermodynamic level (k+1) is between momentum level (k+1)
! and momentum level (k).
tkp1 = level + 1
! Thermodynamic level (k) is between momentum level (k)
! and momentum level (k-1).
tk = level
if ( level == 1 ) then
! k = 1; lower boundery level at surface.
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs(kp1_mdiag) &
= zero
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(k_mdiag) &
= zero
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(km1_mdiag) &
= zero
elseif ( level > 1 .and. level < gr%nz ) then
! Most of the interior model; normal conditions.
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs(kp1_mdiag) &
= + wm_zm * invrs_dzm * gr%weights_zm2zt(m_above,tkp1)
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(k_mdiag) &
= + wm_zm * invrs_dzm * ( gr%weights_zm2zt(m_below,tkp1) &
- gr%weights_zm2zt(m_above,tk) )
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(km1_mdiag) &
= - wm_zm * invrs_dzm * gr%weights_zm2zt(m_below,tk)
elseif ( level == gr%nz ) then
! k = gr%nz (top level); upper boundary level.
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs(kp1_mdiag) &
= zero
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(k_mdiag) &
= zero
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(km1_mdiag) &
= zero
endif
return
end function term_ma_zm_lhs
!=============================================================================================
pure subroutine term_ma_zm_lhs_all( wm_zm, invrs_dzm, & ! Intent(in)
lhs_ma ) ! Intent(out)
! Description:
! This subroutine is an optimized version of term_ma_zm_lhs. term_ma_zm_lhs
! returns a single 3 dimensional array for any specified grid level. This subroutine returns
! an array of 3 dimensional arrays, one for every grid level not including boundary values.
!
! Notes:
! This subroutine exists for performance concerns. It returns all lhs arrays at once
! so that it can be properly vectorized, see clubb:ticket:834 for detail.
!
! THIS SUBROUTINE DOES NOT HANDLE BOUNDARY CONDITIONS AND SETS THEM TO 0
!---------------------------------------------------------------------------------------------
use grid_class, only: &
gr ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
!------------------- Input Variables -------------------
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
wm_zm, & ! wm_zm(k) [m/s]
invrs_dzm ! Inverse of grid spacing (k) [1/m]
!------------------- Output Variables -------------------
real( kind = core_rknd ), dimension(3,gr%nz), intent(out) :: &
lhs_ma
!---------------- Local Variables -------------------
integer :: &
k ! Loop variable for current grid level
!---------------- Begin Code -------------------
! Set lower boundary array to 0
lhs_ma(:,1) = 0.0_core_rknd
! Most of the interior model; normal conditions.
do k = 2, gr%nz-1
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs_ma(1,k) = + wm_zm(k) * invrs_dzm(k) * gr%weights_zm2zt(1,k+1)
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs_ma(2,k) = + wm_zm(k) * invrs_dzm(k) * ( gr%weights_zm2zt(2,k+1) &
- gr%weights_zm2zt(1,k) )
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs_ma(3,k) = - wm_zm(k) * invrs_dzm(k) * gr%weights_zm2zt(2,k)
end do
! Set upper boundary array to 0
lhs_ma(:,gr%nz) = 0.0_core_rknd
return
end subroutine term_ma_zm_lhs_all
end module mean_adv