-
Notifications
You must be signed in to change notification settings - Fork 0
/
advance_helper_module.F90
2007 lines (1533 loc) · 70.4 KB
/
advance_helper_module.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
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module advance_helper_module
! Description:
! This module contains helper methods for the advance_* modules.
!------------------------------------------------------------------------
implicit none
public :: &
set_boundary_conditions_lhs, &
set_boundary_conditions_rhs, &
calc_stability_correction, &
calc_brunt_vaisala_freq_sqd, &
compute_Cx_fnc_Richardson, &
wp2_term_splat_lhs, &
wp3_term_splat_lhs, &
smooth_min, smooth_max, &
smooth_heaviside_peskin, &
calc_xpwp, &
vertical_avg, &
vertical_integral, &
Lscale_width_vert_avg
interface calc_xpwp
module procedure calc_xpwp_1D
module procedure calc_xpwp_2D
end interface
private ! Set Default Scope
!===============================================================================
interface smooth_min
! These functions smooth the output of the min function
! by introducing a varyingly steep path between the two input variables.
! The degree to which smoothing is applied depends on the value of 'smth_coef'.
! If 'smth_coef' goes toward 0, the output of the min function will be
! 0.5 * ((a+b) - abs(a-b))
! If a > b, then this comes out to be b. Likewise, if a < b, abs(a-b)=b-a so we get a.
! Increasing the smoothing coefficient will lead to a greater degree of smoothing
! in the smooth min and max functions. Generally, the coefficient should roughly scale
! with the magnitude of data in the data structure that is to be smoothed, in order to
! obtain a sensible degree of smoothing (not too much, not too little).
module procedure smooth_min_scalar_array
module procedure smooth_min_array_scalar
module procedure smooth_min_arrays
module procedure smooth_min_scalars
end interface
!===============================================================================
interface smooth_max
! These functions smooth the output of the max functions
! by introducing a varyingly steep path between the two input variables.
! The degree to which smoothing is applied depends on the value of 'smth_coef'.
! If 'smth_coef' goes toward 0, the output of the max function will be
! 0.5 * ((a+b) + abs(a-b))
! If a > b, then this comes out to be a. Likewise, if a < b, abs(a-b)=b-a so we get b.
! Increasing the smoothing coefficient will lead to a greater degree of smoothing
! in the smooth min and max functions. Generally, the coefficient should roughly scale
! with the magnitude of data in the data structure that is to be smoothed, in order to
! obtain a sensible degree of smoothing (not too much, not too little).
module procedure smooth_max_scalar_array
module procedure smooth_max_array_scalar
module procedure smooth_max_arrays
module procedure smooth_max_scalars
end interface
!===============================================================================
contains
!---------------------------------------------------------------------------
subroutine set_boundary_conditions_lhs( diag_index, low_bound, high_bound, &
lhs, &
diag_index2, low_bound2, high_bound2 )
! Description:
! Sets the boundary conditions for a left-hand side LAPACK matrix.
!
! References:
! none
!---------------------------------------------------------------------------
use clubb_precision, only: &
core_rknd ! Variable(s)
use constants_clubb, only: &
one, zero
implicit none
! Exernal
intrinsic :: present
! Input Variables
integer, intent(in) :: &
diag_index, low_bound, high_bound ! boundary indexes for the first variable
! Input / Output Variables
real( kind = core_rknd ), dimension(:,:), intent(inout) :: &
lhs ! left hand side of the LAPACK matrix equation
! Optional Input Variables
integer, intent(in), optional :: &
diag_index2, low_bound2, high_bound2 ! boundary indexes for the second variable
! --------------------- BEGIN CODE ----------------------
if ( ( present( low_bound2 ) .or. present( high_bound2 ) ) .and. &
( .not. present( diag_index2 ) ) ) then
error stop "Boundary index provided without diag_index."
end if
! Set the lower boundaries for the first variable
lhs(:,low_bound) = zero
lhs(diag_index,low_bound) = one
! Set the upper boundaries for the first variable
lhs(:,high_bound) = zero
lhs(diag_index,high_bound) = one
! Set the lower boundaries for the second variable, if it is provided
if ( present( low_bound2 ) ) then
lhs(:,low_bound2) = zero
lhs(diag_index2,low_bound2) = one
end if
! Set the upper boundaries for the second variable, if it is provided
if ( present( high_bound2 ) ) then
lhs(:,high_bound2) = zero
lhs(diag_index2,high_bound2) = one
end if
return
end subroutine set_boundary_conditions_lhs
!--------------------------------------------------------------------------
subroutine set_boundary_conditions_rhs( &
low_value, low_bound, high_value, high_bound, &
rhs, &
low_value2, low_bound2, high_value2, high_bound2 )
! Description:
! Sets the boundary conditions for a right-hand side LAPACK vector.
!
! References:
! none
!---------------------------------------------------------------------------
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Exernal
intrinsic :: present
! Input Variables
! The values for the first variable
real( kind = core_rknd ), intent(in) :: low_value, high_value
! The bounds for the first variable
integer, intent(in) :: low_bound, high_bound
! Input / Output Variables
! The right-hand side vector
real( kind = core_rknd ), dimension(:), intent(inout) :: rhs
! Optional Input Variables
! The values for the second variable
real( kind = core_rknd ), intent(in), optional :: low_value2, high_value2
! The bounds for the second variable
integer, intent(in), optional :: low_bound2, high_bound2
! -------------------- BEGIN CODE ------------------------
! Stop execution if a boundary was provided without a value
if ( (present( low_bound2 ) .and. (.not. present( low_value2 ))) .or. &
(present( high_bound2 ) .and. (.not. present( high_value2 ))) ) then
error stop "Boundary condition provided without value."
end if
! Set the lower and upper bounds for the first variable
rhs(low_bound) = low_value
rhs(high_bound) = high_value
! If a lower bound was given for the second variable, set it
if ( present( low_bound2 ) ) then
rhs(low_bound2) = low_value2
end if
! If an upper bound was given for the second variable, set it
if ( present( high_bound2 ) ) then
rhs(high_bound2) = high_value2
end if
return
end subroutine set_boundary_conditions_rhs
!===============================================================================
subroutine calc_stability_correction( nz, ngrdcol, gr, &
thlm, Lscale, em, &
exner, rtm, rcm, &
p_in_Pa, thvm, ice_supersat_frac, &
lambda0_stability_coef, &
bv_efold, &
l_brunt_vaisala_freq_moist, &
l_use_thvm_in_bv_freq, &
stability_correction )
!
! Description:
! Stability Factor
!
! References:
!
!--------------------------------------------------------------------
use constants_clubb, only: &
zero, one, three ! Constant(s)
use grid_class, only: &
grid, & ! Type
zt2zm ! Procedure(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! ---------------- Input Variables ----------------
integer, intent(in) :: &
nz, &
ngrdcol
type (grid), target, intent(in) :: gr
real( kind = core_rknd ), intent(in), dimension(ngrdcol,nz) :: &
Lscale, & ! Turbulent mixing length [m]
em, & ! Turbulent Kinetic Energy (TKE) [m^2/s^2]
thlm, & ! th_l (thermo. levels) [K]
exner, & ! Exner function [-]
rtm, & ! total water mixing ratio, r_t [kg/kg]
rcm, & ! cloud water mixing ratio, r_c [kg/kg]
p_in_Pa, & ! Air pressure [Pa]
thvm, & ! Virtual potential temperature [K]
ice_supersat_frac
real( kind = core_rknd ), intent(in) :: &
lambda0_stability_coef, & ! CLUBB tunable parameter lambda0_stability_coef
bv_efold ! Control parameter for inverse e-folding of
! cloud fraction in the mixed Brunt Vaisala frequency
logical, intent(in) :: &
l_brunt_vaisala_freq_moist, & ! Use a different formula for the Brunt-Vaisala frequency in
! saturated atmospheres (from Durran and Klemp, 1982)
l_use_thvm_in_bv_freq ! Use thvm in the calculation of Brunt-Vaisala frequency
! ---------------- Output Variables ----------------
real( kind = core_rknd ), intent(out), dimension(ngrdcol,nz) :: &
stability_correction
! ---------------- Local Variables ----------------
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
brunt_vaisala_freq_sqd, & ! []
brunt_vaisala_freq_sqd_mixed, &
brunt_vaisala_freq_sqd_dry, & ! []
brunt_vaisala_freq_sqd_moist, &
lambda0_stability, &
Lscale_zm
integer :: i, k
!------------ Begin Code --------------
!$acc enter data create( brunt_vaisala_freq_sqd, brunt_vaisala_freq_sqd_mixed, &
!$acc brunt_vaisala_freq_sqd_moist, brunt_vaisala_freq_sqd_dry, &
!$acc lambda0_stability, Lscale_zm )
call calc_brunt_vaisala_freq_sqd( nz, ngrdcol, gr, thlm, & ! intent(in)
exner, rtm, rcm, p_in_Pa, thvm, & ! intent(in)
ice_supersat_frac, & ! intent(in)
l_brunt_vaisala_freq_moist, & ! intent(in)
l_use_thvm_in_bv_freq, & ! intent(in)
bv_efold, & ! intent(in)
brunt_vaisala_freq_sqd, & ! intent(out)
brunt_vaisala_freq_sqd_mixed,& ! intent(out)
brunt_vaisala_freq_sqd_dry, & ! intent(out)
brunt_vaisala_freq_sqd_moist ) ! intent(out)
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
if ( brunt_vaisala_freq_sqd(i,k) > zero ) then
lambda0_stability(i,k) = lambda0_stability_coef
else
lambda0_stability(i,k) = zero
end if
end do
end do
!$acc end parallel loop
Lscale_zm = zt2zm( nz, ngrdcol, gr, Lscale(:,:) )
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
stability_correction(i,k) = one + min( lambda0_stability(i,k) * brunt_vaisala_freq_sqd(i,k) &
* Lscale_zm(i,k)**2 / em(i,k), three )
end do
end do
!$acc end parallel loop
!$acc exit data delete( brunt_vaisala_freq_sqd, brunt_vaisala_freq_sqd_mixed, &
!$acc brunt_vaisala_freq_sqd_moist, brunt_vaisala_freq_sqd_dry, &
!$acc lambda0_stability, Lscale_zm )
return
end subroutine calc_stability_correction
!===============================================================================
subroutine calc_brunt_vaisala_freq_sqd( nz, ngrdcol, gr, thlm, &
exner, rtm, rcm, p_in_Pa, thvm, &
ice_supersat_frac, &
l_brunt_vaisala_freq_moist, &
l_use_thvm_in_bv_freq, &
bv_efold, &
brunt_vaisala_freq_sqd, &
brunt_vaisala_freq_sqd_mixed,&
brunt_vaisala_freq_sqd_dry, &
brunt_vaisala_freq_sqd_moist )
! Description:
! Calculate the Brunt-Vaisala frequency squared, N^2.
! References:
! ?
!-----------------------------------------------------------------------
use clubb_precision, only: &
core_rknd ! Konstant
use constants_clubb, only: &
grav, & ! Constant(s)
Lv, &
Cp, &
Rd, &
ep, &
one, &
one_half, &
zero_threshold
use parameters_model, only: &
T0 ! Variable!
use grid_class, only: &
grid, & ! Type
ddzt, & ! Procedure(s)
zt2zm
use T_in_K_module, only: &
thlm2T_in_K ! Procedure
use saturation, only: &
sat_mixrat_liq ! Procedure
implicit none
!---------------------------- Input Variables ----------------------------
integer, intent(in) :: &
nz, &
ngrdcol
type (grid), target, intent(in) :: gr
real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
thlm, & ! th_l (thermo. levels) [K]
exner, & ! Exner function [-]
rtm, & ! total water mixing ratio, r_t [kg/kg]
rcm, & ! cloud water mixing ratio, r_c [kg/kg]
p_in_Pa, & ! Air pressure [Pa]
thvm, & ! Virtual potential temperature [K]
ice_supersat_frac
logical, intent(in) :: &
l_brunt_vaisala_freq_moist, & ! Use a different formula for the Brunt-Vaisala frequency in
! saturated atmospheres (from Durran and Klemp, 1982)
l_use_thvm_in_bv_freq ! Use thvm in the calculation of Brunt-Vaisala frequency
real( kind = core_rknd ), intent(in) :: &
bv_efold ! Control parameter for inverse e-folding of
! cloud fraction in the mixed Brunt Vaisala frequency
!---------------------------- Output Variables ----------------------------
real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
brunt_vaisala_freq_sqd, & ! Brunt-Vaisala frequency squared, N^2 [1/s^2]
brunt_vaisala_freq_sqd_mixed, &
brunt_vaisala_freq_sqd_dry,&
brunt_vaisala_freq_sqd_moist
!---------------------------- Local Variables ----------------------------
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
T_in_K, T_in_K_zm, rsat, rsat_zm, thm, thm_zm, ddzt_thlm, &
ddzt_thm, ddzt_rsat, ddzt_rtm, thvm_zm, ddzt_thvm
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
stat_dry, stat_liq, ddzt_stat_liq, ddzt_stat_liq_zm, &
stat_dry_virtual, stat_dry_virtual_zm, ddzt_rtm_zm
integer :: i, k
!---------------------------- Begin Code ----------------------------
!$acc data copyin( gr, gr%zt, &
!$acc thlm, exner, rtm, rcm, p_in_Pa, thvm, ice_supersat_frac ) &
!$acc copyout( brunt_vaisala_freq_sqd, brunt_vaisala_freq_sqd_mixed, &
!$acc brunt_vaisala_freq_sqd_dry, brunt_vaisala_freq_sqd_moist ) &
!$acc create( T_in_K, T_in_K_zm, rsat, rsat_zm, thm, thm_zm, ddzt_thlm, &
!$acc ddzt_thm, ddzt_rsat, ddzt_rtm, thvm_zm, ddzt_thvm, stat_dry, &
!$acc stat_liq, ddzt_stat_liq, ddzt_stat_liq_zm, stat_dry_virtual, &
!$acc stat_dry_virtual_zm, ddzt_rtm_zm )
ddzt_thlm = ddzt( nz, ngrdcol, gr, thlm )
thvm_zm = zt2zm( nz, ngrdcol, gr, thvm )
ddzt_thvm = ddzt( nz, ngrdcol, gr, thvm )
! Dry Brunt-Vaisala frequency
if ( l_use_thvm_in_bv_freq ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
brunt_vaisala_freq_sqd(i,k) = ( grav / thvm_zm(i,k) ) * ddzt_thvm(i,k)
end do
end do
!$acc end parallel loop
else
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
brunt_vaisala_freq_sqd(i,k) = ( grav / T0 ) * ddzt_thlm(i,k)
end do
end do
!$acc end parallel loop
end if
T_in_K = thlm2T_in_K( nz, ngrdcol, thlm, exner, rcm )
T_in_K_zm = zt2zm( nz, ngrdcol, gr, T_in_K )
rsat = sat_mixrat_liq( nz, ngrdcol, p_in_Pa, T_in_K )
rsat_zm = zt2zm( nz, ngrdcol, gr, rsat )
ddzt_rsat = ddzt( nz, ngrdcol, gr, rsat )
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
thm(i,k) = thlm(i,k) + Lv/(Cp*exner(i,k)) * rcm(i,k)
end do
end do
!$acc end parallel loop
thm_zm = zt2zm( nz, ngrdcol, gr, thm )
ddzt_thm = ddzt( nz, ngrdcol, gr, thm )
ddzt_rtm = ddzt( nz, ngrdcol, gr, rtm )
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
stat_dry(i,k) = Cp * T_in_K(i,k) + grav * gr%zt(i,k)
stat_liq(i,k) = stat_dry(i,k) - Lv * rcm(i,k)
end do
end do
!$acc end parallel loop
ddzt_stat_liq = ddzt( nz, ngrdcol, gr, stat_liq )
ddzt_stat_liq_zm = zt2zm( nz, ngrdcol, gr, ddzt_stat_liq)
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
stat_dry_virtual(i,k) = stat_dry(i,k) + Cp * T_in_K(i,k) &
*( 0.608 * ( rtm(i,k) - rcm(i,k) )- rcm(i,k) )
end do
end do
!$acc end parallel loop
stat_dry_virtual_zm = zt2zm( nz, ngrdcol, gr, stat_dry_virtual)
ddzt_rtm_zm = zt2zm( nz, ngrdcol, gr, ddzt_rtm )
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
brunt_vaisala_freq_sqd_dry(i,k) = ( grav / thm_zm(i,k) )* ddzt_thm(i,k)
end do
end do
!$acc end parallel loop
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
! In-cloud Brunt-Vaisala frequency. This is Eq. (36) of Durran and
! Klemp (1982)
brunt_vaisala_freq_sqd_moist(i,k) = &
grav * ( ((one + Lv*rsat_zm(i,k) / (Rd*T_in_K_zm(i,k))) / &
(one + ep*(Lv**2)*rsat_zm(i,k)/(Cp*Rd*T_in_K_zm(i,k)**2))) * &
( (one/thm_zm(i,k) * ddzt_thm(i,k)) + (Lv/(Cp*T_in_K_zm(i,k)))*ddzt_rsat(i,k)) - &
ddzt_rtm(i,k) )
end do
end do ! k=1, gr%nz
!$acc end parallel loop
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
brunt_vaisala_freq_sqd_mixed(i,k) = &
brunt_vaisala_freq_sqd_moist(i,k) + &
exp( - bv_efold * ice_supersat_frac(i,k) ) * &
( brunt_vaisala_freq_sqd_dry(i,k) - brunt_vaisala_freq_sqd_moist(i,k) )
end do
end do
!$acc end parallel loop
if ( l_brunt_vaisala_freq_moist ) then
brunt_vaisala_freq_sqd = brunt_vaisala_freq_sqd_moist
end if
!$acc end data
return
end subroutine calc_brunt_vaisala_freq_sqd
!===============================================================================
subroutine compute_Cx_fnc_Richardson( nz, ngrdcol, gr, &
thlm, um, vm, em, Lscale, exner, rtm, &
rcm, p_in_Pa, thvm, rho_ds_zm, &
ice_supersat_frac, &
clubb_params, &
l_brunt_vaisala_freq_moist, &
l_use_thvm_in_bv_freq, &
l_use_shear_Richardson, &
l_modify_limiters_for_cnvg_test, &
stats_metadata, &
stats_zm, &
Cx_fnc_Richardson )
! Description:
! Compute Cx as a function of the Richardson number
! References:
! cam:ticket:59
!-----------------------------------------------------------------------
use clubb_precision, only: &
core_rknd ! Konstant
use grid_class, only: &
grid, & ! Type
ddzt, & ! Procedure(s)
zt2zm, &
zm2zt2zm
use constants_clubb, only: &
one, zero
use interpolation, only: &
linear_interp_factor ! Procedure
use parameter_indices, only: &
nparams, & ! Variable(s)
iCx_min, &
iCx_max, &
iRichardson_num_min, &
iRichardson_num_max, &
ibv_efold
use stats_variables, only: &
stats_metadata_type
use stats_type_utilities, only: &
stat_update_var ! Procedure
use stats_type, only: stats ! Type
implicit none
!------------------------------ Constant Parameters ------------------------------
real( kind = core_rknd ), parameter :: &
Richardson_num_divisor_threshold = 1.0e-6_core_rknd, &
Cx_fnc_Richardson_below_ground_value = one
logical, parameter :: &
l_Cx_fnc_Richardson_vert_avg = .false. ! Vertically average Cx_fnc_Richardson over a
! distance of Lscale
real( kind = core_rknd ), parameter :: &
min_max_smth_mag = 1.0e-9_core_rknd ! "base" smoothing magnitude before scaling
! for the respective data structure. See
! https://github.com/larson-group/clubb/issues/965#issuecomment-1119816722
! for a plot on how output behaves with varying min_max_smth_mag
!------------------------------ Input Variables ------------------------------
integer, intent(in) :: &
nz, &
ngrdcol
type (grid), target, intent(in) :: gr
real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
thlm, & ! th_l (liquid water potential temperature) [K]
um, & ! u mean wind component (thermodynamic levels) [m/s]
vm, & ! v mean wind component (thermodynamic levels) [m/s]
em, & ! Turbulent Kinetic Energy (TKE) [m^2/s^2]
Lscale, & ! Turbulent mixing length [m]
exner, & ! Exner function [-]
rtm, & ! total water mixing ratio, r_t [kg/kg]
rcm, & ! cloud water mixing ratio, r_c [kg/kg]
p_in_Pa, & ! Air pressure [Pa]
thvm, & ! Virtual potential temperature [K]
rho_ds_zm, & ! Dry static density on momentum levels [kg/m^3]
ice_supersat_frac ! ice cloud fraction
real( kind = core_rknd ), dimension(nparams), intent(in) :: &
clubb_params ! Array of CLUBB's tunable parameters [units vary]
logical, intent(in) :: &
l_brunt_vaisala_freq_moist, & ! Use a different formula for the Brunt-Vaisala frequency in
! saturated atmospheres (from Durran and Klemp, 1982)
l_use_thvm_in_bv_freq, & ! Use thvm in the calculation of Brunt-Vaisala frequency
l_use_shear_Richardson ! Use shear in the calculation of Richardson number
! Flag to activate modifications on limiters for convergence test
! (smoothed max and min for Cx_fnc_Richardson in advance_helper_module.F90)
! (remove the clippings on brunt_vaisala_freq_sqd_smth in mixing_length.F90)
! (reduce threshold on limiters for sqrt_Ri_zm in mixing_length.F90)
logical, intent(in) :: &
l_modify_limiters_for_cnvg_test
type (stats_metadata_type), intent(in) :: &
stats_metadata
!------------------------------ InOut Variable ------------------------------
type (stats), target, dimension(ngrdcol), intent(inout) :: &
stats_zm
!------------------------------ Output Variable ------------------------------
real( kind = core_rknd), dimension(ngrdcol,nz), intent(out) :: &
Cx_fnc_Richardson
!------------------------------ Local Variables ------------------------------
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
brunt_vaisala_freq_sqd, &
brunt_vaisala_freq_sqd_mixed,&
brunt_vaisala_freq_sqd_dry, &
brunt_vaisala_freq_sqd_moist, &
fnc_Richardson, &
fnc_Richardson_clipped, &
fnc_Richardson_smooth, &
Ri_zm, &
ddzt_um, &
ddzt_vm, &
shear_sqd, &
Lscale_zm, &
Cx_fnc_interp, &
Cx_fnc_Richardson_avg
real ( kind = core_rknd ) :: &
invrs_min_max_diff, &
invrs_num_div_thresh
real( kind = core_rknd ) :: &
Richardson_num_max, & ! CLUBB tunable parameter Richardson_num_max
Richardson_num_min, & ! CLUBB tunable parameter Richardson_num_min
Cx_max, & ! CLUBB tunable parameter max of Cx_fnc_Richardson
Cx_min ! CLUBB tunable parameter min of Cx_fnc_Richardson
integer :: smth_type = 1
integer :: i, k
!------------------------------ Begin Code ------------------------------
!$acc enter data create( brunt_vaisala_freq_sqd, brunt_vaisala_freq_sqd_mixed, &
!$acc brunt_vaisala_freq_sqd_dry, brunt_vaisala_freq_sqd_moist, &
!$acc Cx_fnc_interp, &
!$acc Ri_zm, ddzt_um, ddzt_vm, shear_sqd, Lscale_zm, &
!$acc Cx_fnc_Richardson_avg, fnc_Richardson, &
!$acc fnc_Richardson_clipped, fnc_Richardson_smooth )
call calc_brunt_vaisala_freq_sqd( nz, ngrdcol, gr, thlm, & ! intent(in)
exner, rtm, rcm, p_in_Pa, thvm, & ! intent(in)
ice_supersat_frac, & ! intent(in)
l_brunt_vaisala_freq_moist, & ! intent(in)
l_use_thvm_in_bv_freq, & ! intent(in)
clubb_params(ibv_efold), & ! intent(in)
brunt_vaisala_freq_sqd, & ! intent(out)
brunt_vaisala_freq_sqd_mixed,& ! intent(out)
brunt_vaisala_freq_sqd_dry, & ! intent(out)
brunt_vaisala_freq_sqd_moist ) ! intent(out)
Richardson_num_max = clubb_params(iRichardson_num_max)
Richardson_num_min = clubb_params(iRichardson_num_min)
Cx_max = clubb_params(iCx_max)
Cx_min = clubb_params(iCx_min)
invrs_min_max_diff = one / ( Richardson_num_max - Richardson_num_min )
invrs_num_div_thresh = one / Richardson_num_divisor_threshold
Lscale_zm = zt2zm( nz, ngrdcol, gr, Lscale )
! Calculate shear_sqd
ddzt_um = ddzt( nz, ngrdcol, gr, um )
ddzt_vm = ddzt( nz, ngrdcol, gr, vm )
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
shear_sqd(i,k) = ddzt_um(i,k)**2 + ddzt_vm(i,k)**2
end do
end do
!$acc end parallel loop
if ( stats_metadata%l_stats_samp ) then
!$acc update host(shear_sqd)
do i = 1, ngrdcol
call stat_update_var( stats_metadata%ishear_sqd, shear_sqd(i,:), & ! intent(in)
stats_zm(i) ) ! intent(inout)
end do
end if
if ( l_use_shear_Richardson ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Ri_zm(i,k) = max( 1.0e-7_core_rknd, brunt_vaisala_freq_sqd_mixed(i,k) ) &
/ max( shear_sqd(i,k), 1.0e-7_core_rknd )
end do
end do
!$acc end parallel loop
else
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Ri_zm(i,k) = brunt_vaisala_freq_sqd(i,k) * invrs_num_div_thresh
end do
end do
!$acc end parallel loop
end if
! Cx_fnc_Richardson is interpolated based on the value of Richardson_num
! The min function ensures that Cx does not exceed Cx_max, regardless of the
! value of Richardson_num_max.
if ( l_modify_limiters_for_cnvg_test ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
fnc_Richardson(i,k) = ( Ri_zm(i,k) - Richardson_num_min ) * invrs_min_max_diff
end do
end do
fnc_Richardson_clipped = smooth_min( nz, ngrdcol, one, &
fnc_Richardson, &
min_max_smth_mag )
fnc_Richardson_smooth = smooth_max( nz, ngrdcol, zero, &
fnc_Richardson_clipped, &
min_max_smth_mag )
! use smoothed max amd min to achive smoothed profile and avoid discontinuities
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Cx_fnc_interp(i,k) = fnc_Richardson_smooth(i,k) * ( Cx_max - Cx_min ) + Cx_min
end do
end do
Cx_fnc_Richardson = zm2zt2zm( nz, ngrdcol, gr, Cx_fnc_interp )
else ! default method
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Cx_fnc_Richardson(i,k) = ( max(min(Richardson_num_max, Ri_zm(i,k)), Richardson_num_min) &
- Richardson_num_min ) &
* invrs_min_max_diff * ( Cx_max - Cx_min ) + Cx_min
end do
end do
!$acc end parallel loop
end if
if ( l_Cx_fnc_Richardson_vert_avg ) then
Cx_fnc_Richardson = Lscale_width_vert_avg( nz, ngrdcol, gr, smth_type, &
Cx_fnc_Richardson, Lscale_zm, rho_ds_zm, &
Cx_fnc_Richardson_below_ground_value )
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Cx_fnc_Richardson(i,k) = Cx_fnc_Richardson_avg(i,k)
end do
end do
!$acc end parallel loop
end if
! On some compilers, roundoff error can result in Cx_fnc_Richardson being
! slightly outside the range [0,1]. Thus, it is clipped here.
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Cx_fnc_Richardson(i,k) = max( zero, min( one, Cx_fnc_Richardson(i,k) ) )
end do
end do
!$acc end parallel loop
!$acc exit data delete( brunt_vaisala_freq_sqd, brunt_vaisala_freq_sqd_mixed, &
!$acc brunt_vaisala_freq_sqd_dry, brunt_vaisala_freq_sqd_moist, &
!$acc Cx_fnc_interp, Ri_zm, &
!$acc ddzt_um, ddzt_vm, shear_sqd, Lscale_zm, &
!$acc Cx_fnc_Richardson_avg, fnc_Richardson, &
!$acc fnc_Richardson_clipped, fnc_Richardson_smooth )
return
end subroutine compute_Cx_fnc_Richardson
!----------------------------------------------------------------------
!----------------------------------------------------------------------
function Lscale_width_vert_avg( nz, ngrdcol, gr, smth_type, &
var_profile, Lscale_zm, rho_ds_zm, &
var_below_ground_value )&
result (Lscale_width_vert_avg_output)
! Description:
! Averages a profile with a running mean of width Lscale_zm
! References:
! cam:ticket:59
use clubb_precision, only: &
core_rknd ! Precision
use grid_class, only: &
grid ! Type
use constants_clubb, only: &
zero
implicit none
!-------------------------- Input Variables --------------------------
integer, intent(in) :: &
nz, &
ngrdcol, &
smth_type
type (grid), target, intent(in) :: gr
real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
var_profile, & ! Profile on momentum levels
Lscale_zm, & ! Lscale on momentum levels
rho_ds_zm ! Dry static energy on momentum levels!
real( kind = core_rknd ), intent(in) :: &
var_below_ground_value ! Value to use below ground
!-------------------------- Result Variable --------------------------
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
Lscale_width_vert_avg_output ! Vertically averaged profile (on momentum levels)
!-------------------------- Local Variables --------------------------
integer :: &
k, i, & ! Loop variable
k_avg_lower, &
k_avg_upper, &
k_avg
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
one_half_avg_width, &
numer_terms, &
denom_terms
integer :: &
n_below_ground_levels
real( kind = core_rknd ) :: &
numer_integral, & ! Integral in the numerator (see description)
denom_integral ! Integral in the denominator (see description)
!-------------------------- Begin Code --------------------------
!$acc enter data create( one_half_avg_width, numer_terms, denom_terms )
if ( smth_type == 1 ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
one_half_avg_width(i,k) = max( Lscale_zm(i,k), 500.0_core_rknd )
end do
end do
else if (smth_type == 2 ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
one_half_avg_width(i,k) = 60.0_core_rknd
end do
end do
endif
! Pre calculate numerator and denominator terms
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
numer_terms(i,k) = rho_ds_zm(i,k) * gr%dzm(i,k) * var_profile(i,k)
denom_terms(i,k) = rho_ds_zm(i,k) * gr%dzm(i,k)
end do
end do
! For every grid level
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
!-----------------------------------------------------------------------
! Hunt down all vertical levels with one_half_avg_width(k) of gr%zm(k).
!
! Note: Outdated explanation of version that improves CPU performance
! below. Reworked due to it requiring a k dependency. Now we
! begin looking for k_avg_upper and k_avg_lower starting at
! the kth level.
!
! Outdated but potentially useful note:
! k_avg_upper and k_avg_lower can be saved each loop iteration, this
! reduces iterations beacuse the kth values are likely to be within
! one or two grid levels of the k-1th values. Less searching is required
! by starting the search at the previous values and incrementing or
! decrement as needed.
!-----------------------------------------------------------------------
! Determine if k_avg_upper needs to increment
do k_avg_upper = k, nz-1
if ( gr%zm(i,k_avg_upper+1) - gr%zm(i,k) > one_half_avg_width(i,k) ) then
exit
end if
end do
! Determine if k_avg_lower needs to decrement
do k_avg_lower = k, 2, -1
if ( gr%zm(i,k) - gr%zm(i,k_avg_lower-1) > one_half_avg_width(i,k) ) then
exit
end if
end do
! Compute the number of levels below ground to include.
if ( k_avg_lower > 1 ) then
! k=1, the lowest "real" level, is not included in the average, so no
! below-ground levels should be included.
n_below_ground_levels = 0
numer_integral = zero
denom_integral = zero
else
! The number of below-ground levels included is equal to the distance
! below the lowest level spanned by one_half_avg_width(k)
! divided by the distance between vertical levels below ground; the
! latter is assumed to be the same as the distance between the first and
! second vertical levels.
n_below_ground_levels = int( ( one_half_avg_width(i,k)-(gr%zm(i,k)-gr%zm(i,1)) ) / &
( gr%zm(i,2)-gr%zm(i,1) ) )
numer_integral = n_below_ground_levels * denom_terms(i,1) * var_below_ground_value