forked from ESCOMP/CLUBB_CESM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
diffusion.F90
934 lines (805 loc) · 41.2 KB
/
diffusion.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
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
!-----------------------------------------------------------------------
! $Id$
!===============================================================================
module diffusion
! Description:
! Module diffusion computes the eddy diffusion terms for all of the
! time-tendency (prognostic) equations in the CLUBB parameterization. Most of
! the eddy diffusion terms are solved for completely implicitly, and therefore
! become part of the left-hand side of their respective equations. However,
! wp2 and wp3 have an option to use a Crank-Nicholson eddy diffusion scheme,
! which has both implicit and explicit components.
!
! Function diffusion_zt_lhs handles the eddy diffusion terms for the variables
! located at thermodynamic grid levels. These variables are: wp3 and all
! hydrometeor species. The variables um and vm also use the Crank-Nicholson
! eddy-diffusion scheme for their turbulent advection term.
!
! Function diffusion_zm_lhs handles the eddy diffusion 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 :: diffusion_zt_lhs, &
diffusion_zt_lhs_all, &
diffusion_cloud_frac_zt_lhs, &
diffusion_zm_lhs, &
diffusion_zm_lhs_all
contains
!=============================================================================
pure function diffusion_zt_lhs( K_zm, K_zmm1, nu, &
invrs_dzmm1, invrs_dzm, &
invrs_dzt, level ) &
result( lhs )
! Description:
! Vertical eddy diffusion 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 an eddy diffusion term:
!
! + d [ ( K_zm + nu ) * d(var_zt)/dz ] / dz.
!
! This term is usually solved for completely implicitly, such that:
!
! + d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz.
!
! However, when a Crank-Nicholson scheme is used, the eddy diffusion term
! has both implicit and explicit components, such that:
!
! + (1/2) * d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz
! + (1/2) * d [ ( K_zm + nu ) * d( var_zt(t) )/dz ] / dz;
!
! for which the implicit component is:
!
! + (1/2) * d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz.
!
! Note: When the implicit 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 "-".
!
! Timestep index (t) stands for the index of the current timestep, while
! timestep index (t+1) stands for the index of the next timestep, which is
! being advanced to in solving the d(var_zt)/dt equation.
!
! The implicit portion of this term is discretized as follows:
!
! The values of var_zt are found on the thermodynamic levels, while the
! values of K_zm are found on the momentum levels. The derivatives (d/dz)
! of var_zt are taken over the intermediate momentum levels. At the
! intermediate momentum levels, d(var_zt)/dz is multiplied by ( K_zm + nu ).
! Then, the derivative of the whole mathematical expression is taken over
! the central thermodynamic level, which yields the desired result.
!
! --var_ztp1----------------------------------------------- t(k+1)
!
! ==========d(var_zt)/dz==(K_zm+nu)======================== m(k)
!
! --var_zt-------------------d[(K_zm+nu)*d(var_zt)/dz]/dz-- t(k)
!
! ==========d(var_zt)/dz==(K_zmm1+nu)====================== m(k-1)
!
! --var_ztm1----------------------------------------------- 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) )
! invrs_dzm(k) = 1 / ( zt(k+1) - zt(k) )
! invrs_dzm(k-1) = 1 / ( zt(k) - zt(k-1) )
!
! Note: This function only computes the general implicit form:
! + d [ ( K_zm + nu ) * d( var_zt(t+1) )/dz ] / dz.
! For a Crank-Nicholson scheme, the left-hand side result of this
! function will have to be multiplied by (1/2). For a
! Crank-Nicholson scheme, the right-hand side (explicit) component
! needs to be computed by multiplying the left-hand side results by
! (1/2), reversing the sign on each left-hand side element, and then
! multiplying each element by the appropriate var_zt(t) value from
! the appropriate vertical level.
!
!
! Boundary Conditions:
!
! 1) Zero-flux boundary conditions.
! This function is set up to use zero-flux boundary conditions at both
! the lower boundary level and the upper boundary level. The flux, F,
! is the amount of var_zt flowing normal through the boundary per unit
! time per unit surface area. The derivative of the flux effects the
! time-tendency of var_zt, such that:
!
! d(var_zt)/dt = -dF/dz.
!
! For the 2nd-order eddy-diffusion term, +d[(K_zm+nu)*d(var_zt)/dz]/dz,
! the flux is:
!
! F = -(K_zm+nu)*d(var_zt)/dz.
!
! In order to have zero-flux boundary conditions, the derivative of
! var_zt, d(var_zt)/dz, needs to equal 0 at both the lower boundary and
! the upper boundary.
!
! In order to discretize the lower boundary condition, consider a new
! level outside the model (thermodynamic level 0) just below the lower
! boundary level (thermodynamic level 1). 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 the lower boundary level. Therefore, the value of
! d(var_zt)/dz between the level just outside the model and the lower
! boundary level is 0, satisfying the zero-flux boundary condition. The
! other value for d(var_zt)/dz (between thermodynamic level 2 and
! thermodynamic level 1) is taken over the intermediate momentum level
! (momentum level 1), where it is multiplied by the factor
! ( K_zm(1) + nu ). Then, the derivative of the whole expression is
! taken over the central thermodynamic level.
!
! -var_zt(2)-------------------------------------------- t(2)
!
! ==========d(var_zt)/dz==(K_zm(1)+nu)================== m(1)
!
! -var_zt(1)---------------d[(K_zm+nu)*d(var_zt)/dz]/dz- t(1) Boundary
!
! [d(var_zt)/dz = 0]
!
! -[var_zt(0) = var_zt(1)]-----(level outside model)---- t(0)
!
! The result is dependent only on values of K_zm found at momentum
! level 1 and values of var_zt found at thermodynamic levels 1 and 2.
! Thus, it only affects 2 diagonals on the left-hand side matrix.
!
! The same method can be used to discretize the upper boundary by
! considering a new level outside the model just above the upper boundary
! level.
!
! 2) Fixed-point boundary conditions.
! Many equations in the model use fixed-point boundary conditions rather
! than zero-flux boundary conditions. This means that the value of
! var_zt stays the same over the course of the timestep at the lower
! boundary, as well as at the upper boundary.
!
! In order to discretize the boundary conditions for equations requiring
! fixed-point boundary conditions, either:
! a) in the parent subroutine or function (that calls this function),
! loop over all vertical levels from the second-lowest to the
! second-highest, ignoring the boundary levels. Then set the values
! at the boundary levels in the parent subroutine; or
! b) in the parent subroutine or function, loop over all vertical levels
! and then overwrite the results at the boundary levels.
!
! Either way, at the boundary levels, an array with a value of 1 at the
! main diagonal on the left-hand side and with values of 0 at all other
! diagonals on the left-hand side will preserve the right-hand side value
! at that level, thus satisfying the fixed-point boundary conditions.
!
!
! Conservation Properties:
!
! When zero-flux boundary conditions are used, this technique of
! discretizing the eddy diffusion term leads to conservative differencing.
! When conservative differencing is in place, the column totals for each
! column in the left-hand side matrix (for the eddy diffusion term) should
! be equal to 0. This ensures that the total amount of the quantity var_zt
! over the entire vertical domain is being conserved, meaning that nothing
! is lost due to diffusional effects.
!
! To see that this conservation law is satisfied, compute the eddy diffusion
! of var_zt and integrate vertically. In discretized matrix notation (where
! "i" stands for the matrix column and "j" stands for the matrix row):
!
! 0 = Sum_j Sum_i ( 1/invrs_dzt )_i
! ( invrs_dzt * ((K_zm+nu)*invrs_dzm) )_ij (var_zt)_j.
!
! The left-hand side matrix, ( invrs_dzt * ((K_zm+nu)*invrs_dzm) )_ij, is
! partially written below. The sum over i in the above equation removes
! invrs_dzt everywhere from the matrix below. The sum over j leaves the
! column totals that are desired.
!
! Left-hand side matrix contributions from eddy diffusion term; first four
! vertical levels:
!
! -------------------------------------------------------------------------->
!k=1 | +invrs_dzt(k) -invrs_dzt(k) 0
! | *(K_zm(k)+nu) *(K_zm(k)+nu)
! | *invrs_dzm(k) *invrs_dzm(k)
! |
!k=2 | -invrs_dzt(k) +invrs_dzt(k) -invrs_dzt(k)
! | *(K_zm(k-1)+nu) *[ (K_zm(k)+nu) *(K_zm(k)+nu)
! | *invrs_dzm(k-1) *invrs_dzm(k) *invrs_dzm(k)
! | +(K_zm(k-1)+nu)
! | *invrs_dzm(k-1) ]
! |
!k=3 | 0 -invrs_dzt(k) +invrs_dzt(k)
! | *(K_zm(k-1)+nu) *[ (K_zm(k)+nu)
! | *invrs_dzm(k-1) *invrs_dzm(k)
! | +(K_zm(k-1)+nu)
! | *invrs_dzm(k-1) ]
! |
!k=4 | 0 0 -invrs_dzt(k)
! | *(K_zm(k-1)+nu)
! | *invrs_dzm(k-1)
! \ /
!
! Note: The superdiagonal term from level 3 and both the main diagonal and
! superdiagonal terms from level 4 are not shown on this diagram.
!
! Note: The matrix shown is a tridiagonal matrix. For a band diagonal
! matrix (with 5 diagonals), there would be an extra row between each
! of the rows shown and an extra column between each of the columns
! shown. However, for the purposes of the var_zt eddy diffusion
! term, those extra row and column values are all 0, and the
! conservation properties of the matrix aren't effected.
!
! If fixed-point boundary conditions are used, the matrix entries at
! level 1 (k=1) read: 1 0 0; which means that conservative differencing
! is not in play. The total amount of var_zt over the entire vertical
! domain is not being conserved, as amounts of var_zt may be fluxed out
! through the upper boundary or lower boundary through the effects of
! diffusion.
!
! Brian Griffin. April 26, 2008.
! References:
! None
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! 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.
! Input Variables
real( kind = core_rknd ), intent(in) :: &
K_zm, & ! Coef. of eddy diffusivity at momentum level (k) [m^2/s]
K_zmm1, & ! Coef. of eddy diffusivity at momentum level (k-1) [m^2/s
invrs_dzt, & ! Inverse of grid spacing over thermo. level (k) [1/m]
invrs_dzm, & ! Inverse of grid spacing over momentum level (k) [1/m]
invrs_dzmm1 ! Inverse of grid spacing over momentum level (k-1) [1/m]
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
nu ! Background constant coef. of eddy diffusivity [m^2/s]
integer, intent(in) :: &
level ! Thermodynamic level where calculation occurs. [-]
! Return Variable
real( kind = core_rknd ), dimension(3) :: lhs
if ( level == 1 ) then
! k = 1 (bottom level); lower boundary level.
! Only relevant if zero-flux boundary conditions are used.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) = - invrs_dzt * (K_zm+nu(1)) * invrs_dzm
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) = + invrs_dzt * (K_zm+nu(1)) * invrs_dzm
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) = 0.0_core_rknd
elseif ( level > 1 .and. level < gr%nz ) then
! Most of the interior model; normal conditions.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) = - invrs_dzt * (K_zm+nu(level)) * invrs_dzm
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) = + invrs_dzt * ( (K_zm+nu(level))*invrs_dzm &
+ (K_zmm1+nu(level))*invrs_dzmm1 )
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) = - invrs_dzt * (K_zmm1+nu(level)) * invrs_dzmm1
elseif ( level == gr%nz ) then
! k = gr%nz (top level); upper boundary level.
! Only relevant if zero-flux boundary conditions are used.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) = 0.0_core_rknd
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(k_tdiag) = + invrs_dzt * (K_zmm1+nu(gr%nz)) * invrs_dzmm1
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) = - invrs_dzt * (K_zmm1+nu(gr%nz)) * invrs_dzmm1
endif
end function diffusion_zt_lhs
!====================================================================================
pure subroutine diffusion_zt_lhs_all( K_zm, nu, & ! Intent(in)
invrs_dzm, invrs_dzt, & ! Intent(in)
lhs ) ! Intent(out)
! Description:
! This subroutine is an optimized version of diffusion_zt_lhs. diffusion_zt_lhs
! returns a length 3 array for any specified grid level. This subroutine returns
! an array of length 3 arrays, one for every specified grid level.
!
! Notes:
! This subroutine exists for performance concerns. It returns many 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)
implicit none
!------------------- Input Variables -------------------
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
K_zm, & ! Coef. of eddy diffusivity at momentum level (k) [m^2/s]
invrs_dzt, & ! Inverse of grid spacing over thermo. level (k) [1/m]
invrs_dzm, & ! Inverse of grid spacing over momentum level (k) [1/m]
nu ! Background constant coef. of eddy diffusivity [m^2/s]
!------------------- Output Variables -------------------
real( kind = core_rknd ), dimension(3,gr%nz), intent(out) :: &
lhs ! Eddy diffusion elements of the left hand side matrix
!---------------- Local Variables -------------------
integer :: &
k ! Loop variable for current grid level
!---------------- Begin Code -------------------
! Handle lower boundary level
lhs(1,1) = - invrs_dzt(1) * ( K_zm(1) + nu(1) ) * invrs_dzm(1)
lhs(2,1) = + invrs_dzt(1) * ( K_zm(1) + nu(1) ) * invrs_dzm(1)
lhs(3,1) = 0.0_core_rknd
! Handle every other level
do k = 2, gr%nz-1
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(1,k) = - invrs_dzt(k) * ( K_zm(k) + nu(k) ) * invrs_dzm(k)
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
lhs(2,k) = + invrs_dzt(k) * ( ( K_zm(k) + nu(k) ) * invrs_dzm(k) &
+ ( K_zm(k-1) + nu(k) ) * invrs_dzm(k-1) )
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(3,k) = - invrs_dzt(k) * ( K_zm(k-1) + nu(k) ) * invrs_dzm(k-1)
end do
! Handle upper boundary level
lhs(1,gr%nz) = 0.0_core_rknd
lhs(2,gr%nz) = + invrs_dzt(gr%nz) * ( K_zm(gr%nz-1) + nu(gr%nz) ) * invrs_dzm(gr%nz-1)
lhs(3,gr%nz) = - invrs_dzt(gr%nz) * ( K_zm(gr%nz-1) + nu(gr%nz) ) * invrs_dzm(gr%nz-1)
end subroutine diffusion_zt_lhs_all
!=============================================================================
pure function diffusion_cloud_frac_zt_lhs &
( K_zm, K_zmm1, cloud_frac_zt, cloud_frac_ztm1, &
cloud_frac_ztp1, cloud_frac_zm, &
cloud_frac_zmm1, &
nu, invrs_dzmm1, invrs_dzm, invrs_dzt, level ) &
result( lhs )
! Description:
! This function adds a weight of cloud fraction to the existing diffusion
! function for number concentration variables (e.g. cloud droplet number
! concentration). This code should be considered experimental and may
! contain bugs.
! References:
! This algorithm uses equations derived from Guo, et al. 2009.
!-----------------------------------------------------------------------------
use grid_class, only: &
gr ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! External
intrinsic :: min
! Constant parameters
real( kind = core_rknd ), parameter :: &
cf_ratio = 10._core_rknd ! Maximum cloud-fraction coefficient applied to Kh_zm
integer, parameter :: &
kp1_tdiag = 1, & ! Thermodynamic superdiagonal index.
k_tdiag = 2, & ! Thermodynamic main diagonal index.
km1_tdiag = 3 ! Thermodynamic subdiagonal index.
! Input Variables
real( kind = core_rknd ), intent(in) :: &
K_zm, & ! Coef. of eddy diffusivity at mom. level (k) [m^2/s]
K_zmm1, & ! Coef. of eddy diffusivity at mom. level (k-1) [m^2/s]
cloud_frac_zt, & ! Cloud fraction at the thermo. level (k) [-]
cloud_frac_ztm1, & ! Cloud fraction at the thermo. level (k-1) [-]
cloud_frac_ztp1, & ! Cloud fraction at the thermo. level (k+1) [-]
cloud_frac_zm, & ! Cloud fraction at the momentum level (k) [-]
cloud_frac_zmm1, & ! Cloud fraction at the momentum level (k-1) [-]
invrs_dzt, & ! Inverse of grid spacing over thermo. lev. (k) [1/m]
invrs_dzm, & ! Inverse of grid spacing over mom. level (k) [1/m]
invrs_dzmm1 ! Inverse of grid spacing over mom. level (k-1) [1/m]
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
nu ! Background constant coef. of eddy diffusivity [m^2/s]
integer, intent(in) :: &
level ! Thermodynamic level where calculation occurs. [-]
! Return Variable
real( kind = core_rknd ), dimension(3) :: lhs
! ---- Begin Code ----
if ( level == 1 ) then
! k = 1 (bottom level); lower boundary level.
! Only relevant if zero-flux boundary conditions are used.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
! lhs(kp1_tdiag) = - invrs_dzt &
! * (K_zm+nu) &
! * ( cloud_frac_zm / cloud_frac_ztp1 ) * invrs_dzm
lhs(kp1_tdiag) = - invrs_dzt &
* (K_zm &
* min( cloud_frac_zm / cloud_frac_ztp1, cf_ratio ) &
+ nu(1)) * invrs_dzm
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
! lhs(k_tdiag) = + invrs_dzt &
! * (K_zm+nu) &
! * ( cloud_frac_zm / cloud_frac_ztp1 ) * invrs_dzm
lhs(k_tdiag) = + invrs_dzt &
* (K_zm &
* min( cloud_frac_zm / cloud_frac_ztp1, cf_ratio ) &
+ nu(1)) * invrs_dzm
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
lhs(km1_tdiag) = 0.0_core_rknd
else if ( level > 1 .and. level < gr%nz ) then
! Most of the interior model; normal conditions.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
! lhs(kp1_tdiag) = - invrs_dzt &
! * (K_zm+nu) &
! * ( cloud_frac_zm / cloud_frac_ztp1 ) * invrs_dzm
! lhs(kp1_tdiag) = - invrs_dzt &
! * (K_zm &
! * ( cloud_frac_zm / cloud_frac_ztp1 ) &
! + nu ) * invrs_dzm
lhs(kp1_tdiag) = - invrs_dzt &
* (K_zm &
* min( cloud_frac_zm / cloud_frac_ztp1, cf_ratio ) &
+ nu(level) ) * invrs_dzm
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
! lhs(k_tdiag) = + invrs_dzt &
! * ( ((K_zm+nu)*cloud_frac_zm)*invrs_dzm &
! + ((K_zmm1+nu)*cloud_frac_zmm1)*invrs_dzmm1 ) &
! / cloud_frac_zt
! lhs(k_tdiag) = + invrs_dzt &
! * ( nu*(invrs_dzm+invrs_dzmm1) + &
! ( ((K_zm*cloud_frac_zm)*invrs_dzm +
! (K_zmm1*cloud_frac_zmm1)*invrs_dzmm1)&
! / cloud_frac_zt &
! ) &
! )
lhs(k_tdiag) = + invrs_dzt &
* ( nu(level)*(invrs_dzm+invrs_dzmm1) + &
( K_zm*invrs_dzm* &
min( cloud_frac_zm / cloud_frac_zt, &
cf_ratio ) &
+ K_zmm1*invrs_dzmm1* &
min( cloud_frac_zmm1 / cloud_frac_zt, &
cf_ratio ) &
) &
)
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
! lhs(km1_tdiag) = - invrs_dzt * (K_zmm1+nu) * &
! ( cloud_frac_zmm1 / cloud_frac_ztm1 ) * invrs_dzmm1
lhs(km1_tdiag) = - invrs_dzt &
* (K_zmm1 &
* min( cloud_frac_zmm1 / cloud_frac_ztm1, &
cf_ratio ) &
+ nu(level) ) * invrs_dzmm1
else if ( level == gr%nz ) then
! k = gr%nz (top level); upper boundary level.
! Only relevant if zero-flux boundary conditions are used.
! Thermodynamic superdiagonal: [ x var_zt(k+1,<t+1>) ]
lhs(kp1_tdiag) = 0.0_core_rknd
! Thermodynamic main diagonal: [ x var_zt(k,<t+1>) ]
! lhs(k_tdiag) = + invrs_dzt &
! *(K_zmm1+nu) &
! *( cloud_frac_zmm1 / cloud_frac_ztm1 ) * invrs_dzmm1
lhs(k_tdiag) = + invrs_dzt &
* (K_zmm1 &
* min( cloud_frac_zmm1 / cloud_frac_ztm1, &
cf_ratio ) &
+ nu(gr%nz)) * invrs_dzmm1
! Thermodynamic subdiagonal: [ x var_zt(k-1,<t+1>) ]
! lhs(km1_tdiag) = - invrs_dzt * (K_zmm1+nu) * &
! ( cloud_frac_zmm1 / cloud_frac_ztm1 ) * invrs_dzmm1
lhs(km1_tdiag) = - invrs_dzt &
* (K_zmm1 &
* min( cloud_frac_zmm1 / cloud_frac_ztm1, &
cf_ratio ) &
+ nu(gr%nz)) * invrs_dzmm1
end if
return
end function diffusion_cloud_frac_zt_lhs
!=============================================================================
pure function diffusion_zm_lhs( K_zt, K_ztp1, nu, &
invrs_dztp1, invrs_dzt, &
invrs_dzm, level ) &
result( lhs )
! Description:
! Vertical eddy diffusion 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 an eddy diffusion term:
!
! + d [ ( K_zt + nu ) * d(var_zm)/dz ] / dz.
!
! This term is usually solved for completely implicitly, such that:
!
! + d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz.
!
! However, when a Crank-Nicholson scheme is used, the eddy diffusion term
! has both implicit and explicit components, such that:
!
! + (1/2) * d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz
! + (1/2) * d [ ( K_zt + nu ) * d( var_zm(t) )/dz ] / dz;
!
! for which the implicit component is:
!
! + (1/2) * d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz.
!
! Note: When the implicit 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 "-".
!
! Timestep index (t) stands for the index of the current timestep, while
! timestep index (t+1) stands for the index of the next timestep, which is
! being advanced to in solving the d(var_zm)/dt equation.
!
! The implicit portion of this term is discretized as follows:
!
! The values of var_zm are found on the momentum levels, while the values of
! K_zt are found on the thermodynamic levels. The derivatives (d/dz) of
! var_zm are taken over the intermediate thermodynamic levels. At the
! intermediate thermodynamic levels, d(var_zm)/dz is multiplied by
! ( K_zt + nu ). Then, the derivative of the whole mathematical expression
! is taken over the central momentum level, which yields the desired result.
!
! ==var_zmp1=============================================== m(k+1)
!
! ----------d(var_zm)/dz--(K_ztp1+nu)---------------------- t(k+1)
!
! ==var_zm===================d[(K_zt+nu)*d(var_zm)/dz]/dz== m(k)
!
! ----------d(var_zm)/dz--(K_zt+nu)------------------------ t(k)
!
! ==var_zmm1=============================================== 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) )
! invrs_dzt(k+1) = 1 / ( zm(k+1) - zm(k) )
! invrs_dzt(k) = 1 / ( zm(k) - zm(k-1) )
!
! Note: This function only computes the general implicit form:
! + d [ ( K_zt + nu ) * d( var_zm(t+1) )/dz ] / dz.
! For a Crank-Nicholson scheme, the left-hand side result of this
! function will have to be multiplied by (1/2). For a
! Crank-Nicholson scheme, the right-hand side (explicit) component
! needs to be computed by multiplying the left-hand side results by
! (1/2), reversing the sign on each left-hand side element, and then
! multiplying each element by the appropriate var_zm(t) value from
! the appropriate vertical level.
!
!
! Boundary Conditions:
!
! 1) Zero-flux boundary conditions.
! This function is set up to use zero-flux boundary conditions at both
! the lower boundary level and the upper boundary level. The flux, F,
! is the amount of var_zm flowing normal through the boundary per unit
! time per unit surface area. The derivative of the flux effects the
! time-tendency of var_zm, such that:
!
! d(var_zm)/dt = -dF/dz.
!
! For the 2nd-order eddy-diffusion term, +d[(K_zt+nu)*d(var_zm)/dz]/dz,
! the flux is:
!
! F = -(K_zt+nu)*d(var_zm)/dz.
!
! In order to have zero-flux boundary conditions, the derivative of
! var_zm, d(var_zm)/dz, needs to equal 0 at both the lower boundary and
! the upper boundary.
!
! In order to discretize the lower boundary condition, consider a new
! level outside the model (momentum level 0) just below the lower
! boundary level (momentum level 1). The value of var_zm at the level
! just outside the model is defined to be the same as the value of var_zm
! at the lower boundary level. Therefore, the value of d(var_zm)/dz
! between the level just outside the model and the lower boundary level
! is 0, satisfying the zero-flux boundary condition. The other value for
! d(var_zm)/dz (between momentum level 2 and momentum level 1) is taken
! over the intermediate thermodynamic level (thermodynamic level 2),
! where it is multiplied by the factor ( K_zt(2) + nu ). Then, the
! derivative of the whole expression is taken over the central momentum
! level.
!
! =var_zm(2)============================================ m(2)
!
! ----------d(var_zm)/dz==(K_zt(2)+nu)------------------ t(2)
!
! =var_zm(1)===============d[(K_zt+nu)*d(var_zm)/dz]/dz= m(1) Boundary
!
! ----------[d(var_zm)/dz = 0]-------------------------- t(1)
!
! =[var_zm(0) = var_zm(1)]=====(level outside model)==== m(0)
!
! The result is dependent only on values of K_zt found at thermodynamic
! level 2 and values of var_zm found at momentum levels 1 and 2. Thus,
! it only affects 2 diagonals on the left-hand side matrix.
!
! The same method can be used to discretize the upper boundary by
! considering a new level outside the model just above the upper boundary
! level.
!
! 2) Fixed-point boundary conditions.
! Many equations in the model use fixed-point boundary conditions rather
! than zero-flux boundary conditions. This means that the value of
! var_zm stays the same over the course of the timestep at the lower
! boundary, as well as at the upper boundary.
!
! In order to discretize the boundary conditions for equations requiring
! fixed-point boundary conditions, either:
! a) in the parent subroutine or function (that calls this function),
! loop over all vertical levels from the second-lowest to the
! second-highest, ignoring the boundary levels. Then set the values
! at the boundary levels in the parent subroutine; or
! b) in the parent subroutine or function, loop over all vertical levels
! and then overwrite the results at the boundary levels.
!
! Either way, at the boundary levels, an array with a value of 1 at the
! main diagonal on the left-hand side and with values of 0 at all other
! diagonals on the left-hand side will preserve the right-hand side value
! at that level, thus satisfying the fixed-point boundary conditions.
!
!
! Conservation Properties:
!
! When zero-flux boundary conditions are used, this technique of
! discretizing the eddy diffusion term leads to conservative differencing.
! When conservative differencing is in place, the column totals for each
! column in the left-hand side matrix (for the eddy diffusion term) should
! be equal to 0. This ensures that the total amount of the quantity var_zm
! over the entire vertical domain is being conserved, meaning that nothing
! is lost due to diffusional effects.
!
! To see that this conservation law is satisfied, compute the eddy diffusion
! of var_zm and integrate vertically. In discretized matrix notation (where
! "i" stands for the matrix column and "j" stands for the matrix row):
!
! 0 = Sum_j Sum_i ( 1/invrs_dzm )_i
! ( invrs_dzm * ((K_zt+nu)*invrs_dzt) )_ij (var_zm)_j.
!
! The left-hand side matrix, ( invrs_dzm * ((K_zt+nu)*invrs_dzt) )_ij, is
! partially written below. The sum over i in the above equation removes
! invrs_dzm everywhere from the matrix below. The sum over j leaves the
! column totals that are desired.
!
! Left-hand side matrix contributions from eddy diffusion term; first four
! vertical levels:
!
! ---------------------------------------------------------------------->
!k=1 | +invrs_dzm(k) -invrs_dzm(k) 0
! | *(K_zt(k+1)+nu) *(K_zt(k+1)+nu)
! | *invrs_dzt(k+1) *invrs_dzt(k+1)
! |
!k=2 | -invrs_dzm(k) +invrs_dzm(k) -invrs_dzm(k)
! | *(K_zt(k)+nu) *[ (K_zt(k+1)+nu) *(K_zt(k+1)+nu)
! | *invrs_dzt(k) *invrs_dzt(k+1) *invrs_dzt(k+1)
! | +(K_zt(k)+nu)
! | *invrs_dzt(k) ]
! |
!k=3 | 0 -invrs_dzm(k) +invrs_dzm(k)
! | *(K_zt(k)+nu) *[ (K_zt(k+1)+nu)
! | *invrs_dzt(k) *invrs_dzt(k+1)
! | +(K_zt(k)+nu)
! | *invrs_dzt(k) ]
! |
!k=4 | 0 0 -invrs_dzm(k)
! | *(K_zt(k)+nu)
! | *invrs_dzt(k)
! \ /
!
! Note: The superdiagonal term from level 3 and both the main diagonal and
! superdiagonal terms from level 4 are not shown on this diagram.
!
! Note: The matrix shown is a tridiagonal matrix. For a band diagonal
! matrix (with 5 diagonals), there would be an extra row between each
! of the rows shown and an extra column between each of the columns
! shown. However, for the purposes of the var_zm eddy diffusion
! term, those extra row and column values are all 0, and the
! conservation properties of the matrix aren't effected.
!
! If fixed-point boundary conditions are used, the matrix entries at
! level 1 (k=1) read: 1 0 0; which means that conservative differencing
! is not in play. The total amount of var_zm over the entire vertical
! domain is not being conserved, as amounts of var_zm may be fluxed out
! through the upper boundary or lower boundary through the effects of
! diffusion.
!
! Brian Griffin. April 26, 2008.
! References:
! None
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable(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.
! Input Variables
real( kind = core_rknd ), intent(in) :: &
K_zt, & ! Coef. of eddy diffusivity at thermo. level (k) [m^2/s]
K_ztp1, & ! Coef. of eddy diffusivity at thermo. level (k+1) [m^2/s]
invrs_dzm, & ! Inverse of grid spacing over momentum level (k) [1/m]
invrs_dzt, & ! Inverse of grid spacing over thermo. level (k) [1/m]
invrs_dztp1 ! Inverse of grid spacing over thermo. level (k+1) [1/m]
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
nu ! Background constant coef. of eddy diffusivity [m^2/s]
integer, intent(in) :: &
level ! Momentum level where calculation occurs. [-]
! Return Variable
real( kind = core_rknd ), dimension(3) :: lhs
if ( level == 1 ) then
! k = 1; lower boundary level at surface.
! Only relevant if zero-flux boundary conditions are used.
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs(kp1_mdiag) = - invrs_dzm * (K_ztp1+nu(2)) * invrs_dztp1
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(k_mdiag) = + invrs_dzm * (K_ztp1+nu(2)) * invrs_dztp1
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(km1_mdiag) = 0.0_core_rknd
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) = - invrs_dzm * (K_ztp1+nu(level+1)) * invrs_dztp1
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(k_mdiag) = + invrs_dzm * ( (K_ztp1+nu(level+1))*invrs_dztp1 &
+ (K_zt+nu(level))*invrs_dzt )
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(km1_mdiag) = - invrs_dzm * (K_zt+nu(level)) * invrs_dzt
elseif ( level == gr%nz ) then
! k = gr%nz (top level); upper boundary level.
! Only relevant if zero-flux boundary conditions are used.
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs(kp1_mdiag) = 0.0_core_rknd
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(k_mdiag) = + invrs_dzm * (K_zt+nu(gr%nz)) * invrs_dzt
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(km1_mdiag) = - invrs_dzm * (K_zt+nu(gr%nz)) * invrs_dzt
endif
end function diffusion_zm_lhs
!===================================================================================
pure subroutine diffusion_zm_lhs_all( K_zt, nu, & ! Intent(in)
invrs_dzt, invrs_dzm, & ! Intent(in)
lhs ) ! Intent(out
! Description:
! This subroutine is an optimized version of diffusion_zm_lhs. diffusion_zm_lhs
! returns a length 3 array for any specified grid level. This subroutine returns
! an array of length 3 arrays, one for every specified grid level.
!
! Notes:
! This subroutine exists for performance concerns. It returns many 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)
implicit none
!------------------- Input Variables -------------------
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
K_zt, & ! Coef. of eddy diffusivity at thermo. level (k) [m^2/s]
invrs_dzm, & ! Inverse of grid spacing over momentum level (k) [1/m]
invrs_dzt, & ! Inverse of grid spacing over thermo. level (k) [1/m]
nu ! Background constant coef. of eddy diffusivity [m^2/s]
!------------------- Output Variables -------------------
real( kind = core_rknd ), dimension(3,gr%nz), intent(out) :: &
lhs ! Eddy diffusion elements of the left hand side matrix
!---------------- Local Variables -------------------
integer :: &
k ! Loop variable for current grid level
!---------------- Begin Code -------------------
! Handle lower boundary level
lhs(1,1) = - invrs_dzm(1) * ( K_zt(2) + nu(2) ) * invrs_dzt(2)
lhs(2,1) = + invrs_dzm(1) * ( K_zt(2) + nu(2) ) * invrs_dzt(2)
lhs(3,1) = 0.0_core_rknd
! Handle every other level
do k = 2, gr%nz-1
! Momentum superdiagonal: [ x var_zm(k+1,<t+1>) ]
lhs(1,k) = - invrs_dzm(k) * ( K_zt(k+1) + nu(k+1) ) * invrs_dzt(k+1)
! Momentum main diagonal: [ x var_zm(k,<t+1>) ]
lhs(2,k) = + invrs_dzm(k) * ( ( K_zt(k+1) + nu(k+1) ) * invrs_dzt(k+1) &
+ ( K_zt(k) + nu(k) ) * invrs_dzt(k) )
! Momentum subdiagonal: [ x var_zm(k-1,<t+1>) ]
lhs(3,k) = - invrs_dzm(k) * ( K_zt(k) + nu(k) ) * invrs_dzt(k)
end do
! Handle upper boundary level
lhs(1,gr%nz) = 0.0_core_rknd
lhs(2,gr%nz) = + invrs_dzm(gr%nz) * ( K_zt(gr%nz) + nu(gr%nz) ) * invrs_dzt(gr%nz)
lhs(3,gr%nz) = - invrs_dzm(gr%nz) * ( K_zt(gr%nz) + nu(gr%nz) ) * invrs_dzt(gr%nz)
end subroutine diffusion_zm_lhs_all
end module diffusion