forked from ESCOMP/CLUBB_CESM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
adg1_adg2_3d_luhar_pdf.F90
1392 lines (1103 loc) · 58.2 KB
/
adg1_adg2_3d_luhar_pdf.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
! $Id$
!===============================================================================
module adg1_adg2_3d_luhar_pdf
! Description:
! Contains code to calculate the PDF parameters using the Analytic Double
! Gaussian 1 (ADG1) PDF closure, the Analytic Double Gaussian 2 (ADG2) PDF
! closure, or the 3D Luhar PDF closure. The ADG1 PDF is the original PDF used
! by CLUBB in interactive runs.
! References:
!-------------------------------------------------------------------------
implicit none
public :: ADG1_pdf_driver, & ! Procedure(s)
ADG2_pdf_driver, &
Luhar_3D_pdf_driver, &
ADG1_w_closure, &
calc_Luhar_params, &
close_Luhar_pdf
private :: ADG1_ADG2_responder_params, & ! Procedure(s)
backsolve_Luhar_params, &
max_cubic_root
private
contains
!=============================================================================
subroutine ADG1_pdf_driver( wm, rtm, thlm, wp2, rtp2, thlp2, & ! In
Skw, wprtp, wpthlp, sqrt_wp2, & ! In
sigma_sqd_w, mixt_frac_max_mag, & ! In
sclrm, sclrp2, wpsclrp, l_scalar_calc, & ! In
w_1, w_2, rt_1, rt_2, thl_1, thl_2, & ! Out
varnce_w_1, varnce_w_2, varnce_rt_1, & ! Out
varnce_rt_2, varnce_thl_1, varnce_thl_2, & ! Out
mixt_frac, alpha_rt, alpha_thl, & ! Out
sclr_1, sclr_2, varnce_sclr_1, & ! Out
varnce_sclr_2, alpha_sclr ) ! Out
! Calculates the PDF parameters using the Analytic Double Gaussian 1 (ADG1)
! PDF closure.
! References:
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable type(s)
use constants_clubb, only: &
rt_tol, & ! Constant(s)
thl_tol
use parameters_model, only: &
sclr_dim, & ! Variable(s)
sclr_tol
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Input Variables
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
wm, & ! Mean of w-wind comp. (vert. vel.) [m/s]
rtm, & ! Mean of total water mixing ratio [kg/kg]
thlm, & ! Mean of liquid water potential temp. [K]
wp2, & ! Variance of w (overall) [m^2/s^2]
rtp2, & ! Variance of r_t (overall) [(kg/kg)^2]
thlp2, & ! Variance of th_l (overall) [K^2]
Skw, & ! Skewness of w [-]
wprtp, & ! Covariance of w and r_t [(kg/kg)(m/s)]
wpthlp, & ! Covariance of w and th_l [K(m/s)]
sqrt_wp2, & ! Square root of variance of w [m/s]
sigma_sqd_w ! Width of individual w plumes [-]
real( kind = core_rknd ), intent(in) :: &
mixt_frac_max_mag ! Maximum allowable mag. of mixt_frac [-]
real( kind = core_rknd ), dimension(gr%nz, sclr_dim), intent(in) :: &
sclrm, & ! Mean of passive scalar (overall) [units vary]
sclrp2, & ! Variance of passive scalar (overall) [(units vary)^2]
wpsclrp ! Covariance of w and passive scalar [m/s (units vary)]
logical, intent(in) :: &
l_scalar_calc ! Flag to perform calculations for passive scalars
! Output Variables
real( kind = core_rknd ), dimension(gr%nz), intent(out) :: &
w_1, & ! Mean of w (1st PDF component) [m/s]
w_2, & ! Mean of w (2nd PDF component) [m/s]
rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
thl_1, & ! Mean of th_l (1st PDF component) [K]
thl_2, & ! Mean of th_l (2nd PDF component) [K]
varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
mixt_frac, & ! Mixture fraction (weight of 1st PDF component) [-]
alpha_thl, & ! Factor relating to normalized variance for th_l [-]
alpha_rt ! Factor relating to normalized variance for r_t [-]
real( kind = core_rknd ), dimension(gr%nz, sclr_dim), intent(out) :: &
sclr_1, & ! Mean of passive scalar (1st PDF component) [units vary]
sclr_2, & ! Mean of passive scalar (2nd PDF component) [units vary]
varnce_sclr_1, & ! Variance of pass. sclr (1st PDF comp.) [(units vary)^2]
varnce_sclr_2, & ! Variance of pass. sclr (2nd PDF comp.) [(units vary)^2]
alpha_sclr ! Factor relating to normalized variance for sclr [-]
! Local Variables
real( kind = core_rknd ), dimension(gr%nz) :: &
w_1_n, & ! Normalized mean of w (1st PDF component) [-]
w_2_n ! Normalized mean of w (2nd PDF component) [-]
integer :: i ! Loop index
! Calculate the mixture fraction and the PDF component means and variances
! of w.
call ADG1_w_closure( wm, wp2, Skw, sigma_sqd_w, & ! In
sqrt_wp2, mixt_frac_max_mag, & ! In
w_1, w_2, w_1_n, w_2_n, & ! Out
varnce_w_1, varnce_w_2, mixt_frac ) ! Out
! Calculate the PDF component means and variances of rt.
call ADG1_ADG2_responder_params( rtm, rtp2, wp2, sqrt_wp2, & ! In
wprtp, w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, rt_tol, & ! In
rt_1, rt_2, varnce_rt_1, & ! Out
varnce_rt_2, alpha_rt ) ! Out
! Calculate the PDF component means and variances of thl.
call ADG1_ADG2_responder_params( thlm, thlp2, wp2, sqrt_wp2, & ! In
wpthlp, w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, thl_tol, & ! In
thl_1, thl_2, varnce_thl_1, & ! Out
varnce_thl_2, alpha_thl ) ! Out
! Calculate the PDF component means and variances of passive scalars.
if ( l_scalar_calc ) then
do i = 1, sclr_dim
call ADG1_ADG2_responder_params( sclrm(:,i), sclrp2(:,i), & ! In
wp2, sqrt_wp2, wpsclrp(:,i), & ! In
w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, sclr_tol(i), & ! In
sclr_1(:,i), sclr_2(:,i), & ! Out
varnce_sclr_1(:,i), & ! Out
varnce_sclr_2(:,i), & ! Out
alpha_sclr(:,i) ) ! Out
enddo ! i=1, sclr_dim
endif ! l_scalar_calc
return
end subroutine ADG1_pdf_driver
!=============================================================================
subroutine ADG2_pdf_driver( wm, rtm, thlm, wp2, rtp2, thlp2, & ! In
Skw, wprtp, wpthlp, sqrt_wp2, & ! In
sclrm, sclrp2, wpsclrp, l_scalar_calc, & ! In
w_1, w_2, rt_1, rt_2, thl_1, thl_2, & ! Out
varnce_w_1, varnce_w_2, varnce_rt_1, & ! Out
varnce_rt_2, varnce_thl_1, varnce_thl_2, & ! Out
mixt_frac, alpha_rt, alpha_thl, & ! Out
sigma_sqd_w, sclr_1, sclr_2, & ! Out
varnce_sclr_1, varnce_sclr_2, alpha_sclr ) ! Out
! Calculates the PDF parameters using the Analytic Double Gaussian 2 (ADG2)
! PDF closure.
! References:
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable type(s)
use constants_clubb, only: &
one, & ! Constant(s)
w_tol_sqd, &
rt_tol, &
thl_tol
use parameters_model, only: &
sclr_dim, & ! Variable(s)
sclr_tol
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Input Variables
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
wm, & ! Mean of w-wind component (vertical velocity) [m/s]
rtm, & ! Mean of total water mixing ratio [kg/kg]
thlm, & ! Mean of liquid water potential temperature [K]
wp2, & ! Variance of w (overall) [m^2/s^2]
rtp2, & ! Variance of r_t (overall) [(kg/kg)^2]
thlp2, & ! Variance of th_l (overall) [K^2]
Skw, & ! Skewness of w [-]
wprtp, & ! Covariance of w and r_t [(kg/kg)(m/s)]
wpthlp, & ! Covariance of w and th_l [K(m/s)]
sqrt_wp2 ! Square root of variance of w [m/s]
real( kind = core_rknd ), dimension(gr%nz, sclr_dim), intent(in) :: &
sclrm, & ! Mean of passive scalar (overall) [units vary]
sclrp2, & ! Variance of passive scalar (overall) [(units vary)^2]
wpsclrp ! Covariance of w and passive scalar [m/s (units vary)]
logical, intent(in) :: &
l_scalar_calc ! Flag to perform calculations for passive scalars
! Output Variables
real( kind = core_rknd ), dimension(gr%nz), intent(out) :: &
w_1, & ! Mean of w (1st PDF component) [m/s]
w_2, & ! Mean of w (2nd PDF component) [m/s]
rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
thl_1, & ! Mean of th_l (1st PDF component) [K]
thl_2, & ! Mean of th_l (2nd PDF component) [K]
varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
mixt_frac, & ! Mixture fraction (weight of 1st PDF component) [-]
alpha_thl, & ! Factor relating to normalized variance for th_l [-]
alpha_rt, & ! Factor relating to normalized variance for r_t [-]
sigma_sqd_w ! Width of individual w plumes (ADG1 parameter) [-]
real( kind = core_rknd ), dimension(gr%nz, sclr_dim), intent(out) :: &
sclr_1, & ! Mean of passive scalar (1st PDF component) [units vary]
sclr_2, & ! Mean of passive scalar (2nd PDF component) [units vary]
varnce_sclr_1, & ! Variance of pass. sclr (1st PDF comp.) [(units vary)^2]
varnce_sclr_2, & ! Variance of pass. sclr (2nd PDF comp.) [(units vary)^2]
alpha_sclr ! Factor relating to normalized variance for sclr [-]
! Local Variables
real( kind = core_rknd ), dimension(gr%nz) :: &
w_1_n, & ! Normalized mean of w (1st PDF component) [-]
w_2_n, & ! Normalized mean of w (2nd PDF component) [-]
small_m_w, & ! Luhar's m (tunable parameter) [-]
big_m_w, & ! Luhar's M (a function of Luhar's m) [-]
sigma_sqd_w_1, & ! Normalized width parameter of w (1st PDF component) [-]
sigma_sqd_w_2 ! Normalized width parameter of w (2nd PDF component) [-]
integer :: i ! Loop index
! Calculate the mixture fraction and the PDF component means and variances
! of w.
! Reproduce ADG2_w_closure using separate functions
call calc_Luhar_params( gr%nz, Skw, wp2, & ! intent(in)
wp2, w_tol_sqd, & ! intent(in)
mixt_frac, big_m_w, small_m_w ) ! intent(out)
call close_Luhar_pdf( gr%nz, wm, wp2, mixt_frac, & ! intent(in)
small_m_w, wp2, w_tol_sqd, & ! intent(in)
sigma_sqd_w_1, sigma_sqd_w_2, & ! intent(out)
varnce_w_1, varnce_w_2, & ! intent(out)
w_1_n, w_2_n, w_1, w_2 ) ! intent(out)
! Overwrite sigma_sqd_w for consistency with ADG1
sigma_sqd_w = min( one / ( one + small_m_w**2 ), 0.99_core_rknd )
! Calculate the PDF component means and variances of rt.
call ADG1_ADG2_responder_params( rtm, rtp2, wp2, sqrt_wp2, & ! In
wprtp, w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, rt_tol, & ! In
rt_1, rt_2, varnce_rt_1, & ! Out
varnce_rt_2, alpha_rt ) ! Out
! Calculate the PDF component means and variances of thl.
call ADG1_ADG2_responder_params( thlm, thlp2, wp2, sqrt_wp2, & ! In
wpthlp, w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, thl_tol, & ! In
thl_1, thl_2, varnce_thl_1, & ! Out
varnce_thl_2, alpha_thl ) ! Out
! Calculate the PDF component means and variances of passive scalars.
if ( l_scalar_calc ) then
do i = 1, sclr_dim
call ADG1_ADG2_responder_params( sclrm(:,i), sclrp2(:,i), & ! In
wp2, sqrt_wp2, wpsclrp(:,i), & ! In
w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, sclr_tol(i), & ! In
sclr_1(:,i), sclr_2(:,i), & ! Out
varnce_sclr_1(:,i), & ! Out
varnce_sclr_2(:,i), & ! Out
alpha_sclr(:,i) ) ! Out
enddo ! i=1, sclr_dim
endif ! l_scalar_calc
return
end subroutine ADG2_pdf_driver
!=============================================================================
subroutine Luhar_3D_pdf_driver( wm, rtm, thlm, wp2, rtp2, thlp2, & ! In
Skw, Skrt, Skthl, wprtp, wpthlp, & ! In
w_1, w_2, rt_1, rt_2, thl_1, thl_2, & ! Out
varnce_w_1, varnce_w_2, varnce_rt_1, & ! Out
varnce_rt_2, varnce_thl_1, & ! Out
varnce_thl_2, mixt_frac ) ! Out
! Calculates the PDF parameters using the 3D Luhar PDF closure.
! References:
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable type(s)
use constants_clubb, only: &
w_tol_sqd, & ! Constant(s)
rt_tol, &
thl_tol
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Input Variables
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
wm, & ! Mean of w-wind component (vertical velocity) [m/s]
rtm, & ! Mean of total water mixing ratio [kg/kg]
thlm, & ! Mean of liquid water potential temperature [K]
wp2, & ! Variance of w (overall) [m^2/s^2]
rtp2, & ! Variance of r_t (overall) [(kg/kg)^2]
thlp2, & ! Variance of th_l (overall) [K^2]
Skw, & ! Skewness of w [-]
Skrt, & ! Skewness of rt [-]
Skthl, & ! Skewness of th_l [-]
wprtp, & ! Covariance of w and r_t [(kg/kg)(m/s)]
wpthlp ! Covariance of w and th_l [K(m/s)]
! Output Variables
real( kind = core_rknd ), dimension(gr%nz), intent(out) :: &
w_1, & ! Mean of w (1st PDF component) [m/s]
w_2, & ! Mean of w (2nd PDF component) [m/s]
rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
thl_1, & ! Mean of th_l (1st PDF component) [K]
thl_2, & ! Mean of th_l (2nd PDF component) [K]
varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
mixt_frac ! Mixture fraction (weight of 1st PDF component) [-]
real( kind = core_rknd ), dimension(gr%nz) :: &
w_1_n, & ! Normalized mean of w (1st PDF component) [-]
w_2_n, & ! Normalized mean of w (2nd PDF component) [-]
rt_1_n, & ! Normalized mean of rt (1st PDF component) [-]
rt_2_n, & ! Normalized mean of rt (2nd PDF component) [-]
thl_1_n, & ! Normalized mean of thl (1st PDF component) [-]
thl_2_n, & ! Normalized mean of thl (2nd PDF component) [-]
small_m_w, & ! Luhar's m (tunable parameter) for w [-]
big_m_w, & ! Luhar's M (a function of Luhar's m) for w [-]
small_m_rt, & ! Luhar's m (tunable parameter) for rt [-]
big_m_rt, & ! Luhar's M (a function of Luhar's m) for rt [-]
small_m_thl, & ! Luhar's m (tunable parameter) for thl [-]
big_m_thl, & ! Luhar's M (a function of Luhar's m) for thl [-]
sigma_sqd_w_1, & ! Normalized width parameter of w (1st PDF comp.) [-]
sigma_sqd_w_2, & ! Normalized width parameter of w (2nd PDF comp.) [-]
sigma_sqd_rt_1, & ! Normalized width parameter of rt (1st PDF comp.) [-]
sigma_sqd_rt_2, & ! Normalized width parameter of rt (2nd PDF comp.) [-]
sigma_sqd_thl_1, & ! Normalized width parameter of thl (1st PDF comp.) [-]
sigma_sqd_thl_2 ! Normalized width parameter of thl (2nd PDF comp.) [-]
integer :: k ! Vertical loop index
do k = 1, gr%nz, 1
if ( ( abs( Skw(k) ) >= abs( Skthl(k) ) ) &
.and. ( abs( Skw(k) ) >= abs( Skrt(k) ) ) ) then
! w has the greatest magnitude of skewness.
! Solve for the w PDF
call calc_Luhar_params( 1, Skw(k), wp2(k), & ! In
wp2(k), w_tol_sqd, & ! In
mixt_frac(k), big_m_w(k), small_m_w(k) ) ! Out
call close_Luhar_pdf( 1, wm(k), wp2(k), mixt_frac(k), & ! In
small_m_w(k), wp2(k), w_tol_sqd, & ! In
sigma_sqd_w_1(k), sigma_sqd_w_2(k), & ! Out
varnce_w_1(k), varnce_w_2(k), & ! Out
w_1_n(k), w_2_n(k), w_1(k), w_2(k) ) ! Out
! Solve for the thl PDF
call backsolve_Luhar_params( Skw(k), Skthl(k), & ! In
big_m_w(k), mixt_frac(k), & ! In
big_m_thl(k), small_m_thl(k) ) ! Out
call close_Luhar_pdf( 1, thlm(k), thlp2(k), mixt_frac(k), & ! In
small_m_thl(k), wpthlp(k), thl_tol**2, & ! In
sigma_sqd_thl_1(k), sigma_sqd_thl_2(k), & ! Out
varnce_thl_1(k), varnce_thl_2(k), & ! Out
thl_1_n(k), thl_2_n(k), & ! Out
thl_1(k), thl_2(k) ) ! Out
! Solve for the rt PDF
call backsolve_Luhar_params( Skw(k), Skrt(k), & ! In
big_m_w(k), mixt_frac(k), & ! In
big_m_rt(k), small_m_rt(k) ) ! Out
call close_Luhar_pdf( 1, rtm(k), rtp2(k), mixt_frac(k), & ! In
small_m_rt(k), wprtp(k), rt_tol**2, & ! In
sigma_sqd_rt_1(k), sigma_sqd_rt_2(k), & ! Out
varnce_rt_1(k), varnce_rt_2(k), & ! Out
rt_1_n(k), rt_2_n(k), rt_1(k), rt_2(k) ) ! Out
elseif ( ( abs( Skthl(k) ) > abs( Skw(k) ) ) &
.and. ( abs( Skthl(k) ) >= abs( Skrt(k) ) ) ) then
! theta-l has the greatest magnitude of skewness.
! Solve for the thl PDF
call calc_Luhar_params( 1, Skthl(k), wpthlp(k), & ! In
thlp2(k), thl_tol**2, & ! In
mixt_frac(k), big_m_thl(k), & ! Out
small_m_thl(k) ) ! Out
! Solve for the thl PDF
call close_Luhar_pdf( 1, thlm(k), thlp2(k), mixt_frac(k), & ! In
small_m_thl(k), wpthlp(k), thl_tol**2, & ! In
sigma_sqd_thl_1(k), sigma_sqd_thl_2(k), & ! Out
varnce_thl_1(k), varnce_thl_2(k), & ! Out
thl_1_n(k), thl_2_n(k), & ! Out
thl_1(k), thl_2(k) ) ! Out
! Solve for the w PDF
call backsolve_Luhar_params( Skthl(k), Skw(k), & ! In
big_m_thl(k), mixt_frac(k), & ! In
big_m_w(k), small_m_w(k) ) ! Out
call close_Luhar_pdf( 1, wm(k), wp2(k), mixt_frac(k), & ! In
small_m_w(k), wp2(k), w_tol_sqd, & ! In
sigma_sqd_w_1(k), sigma_sqd_w_2(k), & ! Out
varnce_w_1(k), varnce_w_2(k), & ! Out
w_1_n(k), w_2_n(k), w_1(k), w_2(k) ) ! Out
! Solve for the rt PDF
call backsolve_Luhar_params( Skthl(k), Skrt(k), & ! In
big_m_thl(k), mixt_frac(k), & ! In
big_m_rt(k), small_m_rt(k) ) ! Out
call close_Luhar_pdf( 1, rtm(k), rtp2(k), mixt_frac(k), & ! In
small_m_rt(k), wprtp(k), rt_tol**2, & ! In
sigma_sqd_rt_1(k), sigma_sqd_rt_2(k), & ! Out
varnce_rt_1(k), varnce_rt_2(k), & ! Out
rt_1_n(k), rt_2_n(k), rt_1(k), rt_2(k) ) ! Out
else
! rt has the greatest magnitude of skewness.
! Solve for the rt PDF
call calc_Luhar_params( 1, Skrt(k), wprtp(k), & ! In
rtp2(k), rt_tol**2, & ! In
mixt_frac(k), big_m_rt(k), & ! Out
small_m_rt(k) ) ! Out
! Solve for the rt PDF
call close_Luhar_pdf( 1, rtm(k), rtp2(k), mixt_frac(k), & ! In
small_m_rt(k), wprtp(k), rt_tol**2, & ! In
sigma_sqd_rt_1(k), sigma_sqd_rt_2(k), & ! Out
varnce_rt_1(k), varnce_rt_2(k), & ! Out
rt_1_n(k), rt_2_n(k), rt_1(k), rt_2(k) ) ! Out
! Solve for the w PDF
call backsolve_Luhar_params( Skrt(k), Skw(k), & ! In
big_m_rt(k), mixt_frac(k), & ! In
big_m_w(k), small_m_w(k) ) ! Out
call close_Luhar_pdf( 1, wm(k), wp2(k), mixt_frac(k), & ! In
small_m_w(k), wp2(k), w_tol_sqd, & ! In
sigma_sqd_w_1(k), sigma_sqd_w_2(k), & ! Out
varnce_w_1(k), varnce_w_2(k), & ! Out
w_1_n(k), w_2_n(k), w_1(k), w_2(k) ) ! OUt
! Solve for the thl PDF
call backsolve_Luhar_params( Skrt(k), Skthl(k), & ! In
big_m_rt(k), mixt_frac(k), & ! In
big_m_thl(k), small_m_thl(k) ) ! Out
call close_Luhar_pdf( 1, thlm(k), thlp2(k), mixt_frac(k), & ! In
small_m_thl(k), wpthlp(k), thl_tol**2, & ! In
sigma_sqd_thl_1(k), sigma_sqd_thl_2(k), & ! Out
varnce_thl_1(k), varnce_thl_2(k), & ! Out
thl_1_n(k), thl_2_n(k), & ! Out
thl_1(k), thl_2(k) ) ! Out
endif
enddo ! k = 1, gr%nz, 1
return
end subroutine Luhar_3D_pdf_driver
!=============================================================================
subroutine ADG1_w_closure( wm, wp2, Skw, sigma_sqd_w, & ! In
sqrt_wp2, mixt_frac_max_mag, & ! In
w_1, w_2, w_1_n, w_2_n, & ! Out
varnce_w_1, varnce_w_2, mixt_frac ) ! Out
! Description:
! Calculates the mixture fraction, the PDF component means of w, and the PDF
! component variances of w for the Analytic Double Gaussian 1 (ADG1)
! closure. It sets the widths of both w Gaussians to be the same, and
! furthermore, relates them to the overall variance of w (<w'^2>) by a
! parameter, sigma_sqd_w. The equation is:
!
! sigma_w_1^2 = sigma_w_2^2 = sigma_sqd_w * <w'^2>;
!
! where sigma_w_1^2 is the variance of w in the 1st PDF component,
! sigma_w_2^2 is the variance of w in the 2nd PDF component, and parameter
! sigma_sqd_w must have a value within the range 0 <= sigma_sqd_w <= 1.
!
! References:
! Golaz, J-C., V. E. Larson, and W. R. Cotton, 2002a: A PDF-based model for
! boundary layer clouds. Part I: Method and model description. J. Atmos.
! Sci., 59, 3540–3551.
!
! Vincent E. Larson and Jean-Christophe Golaz, 2005: Using Probability
! Density Functions to Derive Consistent Closure Relationships among
! Higher-Order Moments. Mon. Wea. Rev., 133, 1023–1042.
!-----------------------------------------------------------------------
use grid_class, only: &
gr ! Variable type(s)
use constants_clubb, only: &
four, & ! Constant(s)
one, &
one_half, &
zero, &
w_tol_sqd
use clubb_precision, only: &
core_rknd ! Precision
implicit none
! Input Variables
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
wm, & ! Mean of w (overall) [m/s]
wp2, & ! Variance of w (overall) [m^2/s^2]
Skw, & ! Skewness of w [-]
sigma_sqd_w, & ! Widths of each w Gaussian [-]
sqrt_wp2 ! Square root of the variance of w [m/s]
real( kind = core_rknd ), intent(in) :: &
mixt_frac_max_mag ! Maximum allowable mag. of mixt_frac [-]
! Output Variables
real( kind = core_rknd ), dimension(gr%nz), intent(out) :: &
w_1, & ! Mean of w (1st PDF component) [m/s]
w_2, & ! Mean of w (2nd PDF component) [m/s]
w_1_n, & ! Normalized mean of w (1st PDF component) [-]
w_2_n, & ! Normalized mean of w (2nd PDF component) [-]
varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
mixt_frac ! Mixture fraction [-]
!----- Begin Code -----
where ( wp2 > w_tol_sqd )
! Width (standard deviation) parameters are non-zero
! The variable "mixt_frac" is the weight of the 1st PDF component. The
! weight of the 2nd PDF component is "1-mixt_frac". If there isn't any
! skewness of w (Sk_w = 0 because w'^3 = 0), mixt_frac = 0.5, and both
! PDF components are equally weighted. If there is positive skewness of
! w (Sk_w > 0 because w'^3 > 0), 0 < mixt_frac < 0.5, and the 2nd PDF
! component has greater weight than does the 1st PDF component. If there
! is negative skewness of w (Sk_w < 0 because w'^3 < 0),
! 0.5 < mixt_frac < 1, and the 1st PDF component has greater weight than
! does the 2nd PDF component.
where ( abs( Skw ) <= 1.0e-5_core_rknd )
mixt_frac = one_half
elsewhere
mixt_frac &
= one_half &
* ( one - Skw / sqrt( four * ( one - sigma_sqd_w )**3 + Skw**2 ) )
endwhere
! Clip mixt_frac, and 1 - mixt_frac, to avoid dividing by a small number.
! Formula for mixt_frac_max_mag =
! 1 - ( 1/2 * ( 1 - Skw_max
! / sqrt( 4*( 1 - sigma_sqd_w )^3 + Skw_max^2 ) ) ),
! where sigma_sqd_w is fixed at 0.4 for this calculation.
mixt_frac = min( max( mixt_frac, one - mixt_frac_max_mag ), &
mixt_frac_max_mag )
! The normalized mean of w for Gaussian "plume" 1 is w_1_n. It's value
! will always be greater than 0. As an example, a value of 1.0 would
! indicate that the actual mean of w for Gaussian "plume" 1 is found 1.0
! standard deviation above the overall mean for w.
w_1_n &
= sqrt( ( ( one - mixt_frac ) / mixt_frac ) * ( one - sigma_sqd_w ) )
! The normalized mean of w for Gaussian "plume" 2 is w_2_n. It's value
! will always be less than 0. As an example, a value of -0.5 would
! indicate that the actual mean of w for Gaussian "plume" 2 is found 0.5
! standard deviations below the overall mean for w.
w_2_n &
= -sqrt( ( mixt_frac / ( one - mixt_frac ) ) * ( one - sigma_sqd_w ) )
! The mean of w for Gaussian "plume" 1 is w_1.
w_1 = wm + sqrt_wp2 * w_1_n
! The mean of w for Gaussian "plume" 2 is w_2.
w_2 = wm + sqrt_wp2 * w_2_n
! The variance of w for Gaussian "plume" 1 for varnce_w_1.
varnce_w_1 = sigma_sqd_w * wp2
! The variance of w for Gaussian "plume" 2 for varnce_w_2.
! The variance in both Gaussian "plumes" is defined to be the same.
varnce_w_2 = sigma_sqd_w * wp2
elsewhere
! Vertical velocity doesn't vary.
mixt_frac = one_half
w_1_n = sqrt( one - sigma_sqd_w )
w_2_n = -sqrt( one - sigma_sqd_w )
w_1 = wm
w_2 = wm
varnce_w_1 = zero
varnce_w_2 = zero
endwhere ! Widths non-zero
return
end subroutine ADG1_w_closure
!=============================================================================
subroutine calc_Luhar_params( nz, Skx, wpxp, & ! In
xp2, x_tol_sqd, & ! In
mixt_frac, big_m, small_m ) ! Out
! Description:
! For the Luhar closure, this subroutine takes Skx (and Skw) as input and
! outputs the mixture fraction, big_m, and small_m. This code was written
! using the equations and nomenclature of Larson et al. (2002) Appendix
! section e.
!
! The relationship between skewness of x (Skx), mixture fraction (a), and
! Luhar's small m (m) is given by:
!
! Skx^2 = ( m^2 * ( m^2 + 3 )^2 / ( m^2 + 1 )^3 )
! * ( 1 - 2*a )^2 / ( a * ( 1 - a ) ).
!
! Luhar's large M (M) is used to more easily express the factor involving
! the m's:
!
! M = ( m^2 + 1 )^3 / ( m^2 * ( m^2 + 3 )^2 ).
!
! The equation involving skewness of x becomes:
!
! Skx^2 = ( 1 / M ) * ( 1 - 2*a )^2 / ( a * ( 1 - a ) );
!
! or:
!
! M * Skx^2 = ( 1 - 2*a )^2 / ( a * ( 1 - a ) ).
!
! This equation can be rewritten as:
!
! ( a * ( 1 - a ) ) * M * Skx^2 = ( 1 - 2*a )^2;
!
! as well as:
!
! ( a - a^2 ) * M * Skx^2 = 1 - 4*a + 4*a^2;
!
! and eventually as:
!
! ( 4 + M * Skx^2 ) * a^2 - ( 4 + M * Skx^2 ) * a + 1 = 0.
!
! Solving the quadratic equation for a:
!
! a = (1/2) * ( 1 +- Skx * sqrt( 1 / ( 4/M + Skx^2 ) ) ).
!
! Since by definition, mu_w_1 >= mu_w_2, a < 0.5 when Skw > 0, the equation
! for mixture fraction is:
!
! a = (1/2) * ( 1 - Skx * sqrt( 1 / ( 4/M + Skx^2 ) ) ).
!
! For 3-D Luhar, the variable (w, rt, or theta-l) with the greatest
! magnitude of skewness is used to calculate mixture fraction. Since it is
! desirable to still have a < 0.5 when Skw > 0 and a > 0.5 when Skw < 0, the
! sign function is used. The value of Skx is replaced by:
!
! Skx|_adj = sgn( <w'x'> ) * Skx;
!
! where
!
! sgn( <w'x'> ) = | 1 when <w'x'> >= 0
! | -1 when <w'x'> < 0.
!
! Since Skx|_adj^2 = ( sgn( <w'x'> ) * Skx )^2 = ( sgn( <w'x'> ) )^2 * Skx^2
! = Skx^2, the equation for mixture fraction is:
!
! a = (1/2) * ( 1 - sgn( <w'x'> ) * Skx * sqrt( 1 / ( 4/M + Skx^2 ) ) ).
!
! When using the ADG2 closure or when using the 3-D Luhar closure when the
! variable with the greatest magnitude of skewness is w, Skw = Skx and
! sgn( <w'^2> ) is always equal to 1, reducing the equation to its previous
! form.
! References:
! Vincent E. Larson, Jean-Christophe Golaz, and William R. Cotton, 2002:
! Small-Scale and Mesoscale Variability in Cloudy Boundary Layers: Joint
! Probability Density Functions. J. Atmos. Sci., 59, 3519–3539.
!-----------------------------------------------------------------------
use constants_clubb, only: &
four, & ! Constant(s)
three, &
one, &
two_thirds, &
one_half, &
one_third, &
zero
use clubb_precision, only: &
core_rknd ! Precision
implicit none
! Input Variables
integer, intent(in) :: &
nz ! Number of vertical levels
real( kind = core_rknd ), dimension(nz), intent(in) :: &
Skx, & ! Skewness of x ( <x'^3> / <x'2>^(3/2) ) [-]
wpxp, & ! Covariance of w and x [m/s (x units)]
xp2 ! Variance of x [(x units)^2]
real( kind = core_rknd ), intent(in) :: &
x_tol_sqd ! Tolerance value of x squared [(x units)^2]
! Output Variables
real( kind = core_rknd ), dimension(nz), intent(out) :: &
mixt_frac, & ! Mixture fraction [-]
big_m, & ! Luhar's M [-]
small_m ! Luhar's m [-]
! Local Variables
real( kind = core_rknd ), dimension(nz) :: &
small_m_sqd, & ! Luhar's m^2 [-]
sgn_wpxp ! Sgn(<w'x'>); 1 when <w'x'> >= 0 or -1 when <w'x'> < 0 [-]
where ( xp2 > x_tol_sqd )
! Width (standard deviation) parameters are non-zero
! Calculate Luhar's m (small m).
! If Skx is very small, then small_m will tend to zero which risks
! divide-by-zero. To ameliorate this problem, we enforce abs( x_1_n )
! and abs( x_2_n ) > 0.05.
! Note: Luhar's small_m (m) is the only tunable parameter in the Luhar
! closure, so this equation can be changed. However, the value
! of m should go toward 0 as Skx goes toward 0 so that the double
! Gaussian reduces to a single Gaussian when the distribution is
! unskewed.
small_m = max( two_thirds * abs( Skx )**one_third, 0.05_core_rknd )
! Calculate m^2.
small_m_sqd = small_m**2
! Calculate Luhar's M (big M).
big_m = ( one + small_m_sqd )**3 &
/ ( ( three + small_m_sqd )**2 * small_m_sqd )
! Calculate sgn( <w'x'> ).
where ( wpxp >= zero )
sgn_wpxp = one
elsewhere ! <w'x'> < 0
sgn_wpxp = -one
endwhere ! <w'x'> >= 0
! Calculate mixture fraction.
mixt_frac = one_half &
* ( one - sgn_wpxp * Skx &
* sqrt( one / ( ( four / big_m ) + Skx**2 ) ) )
elsewhere
! Variable x doesn't vary.
mixt_frac = one_half
! For output purposes, set small_m and big_m to 0.
small_m = zero
big_m = zero
endwhere ! Widths non-zero
return
end subroutine calc_Luhar_params
!=============================================================================
subroutine close_Luhar_pdf( nz, xm, xp2, mixt_frac, & ! In
small_m, wpxp, x_tol_sqd, & ! In
sigma_sqd_x_1, sigma_sqd_x_2, & ! Out
varnce_x_1, varnce_x_2, & ! Out
x_1_n, x_2_n, x_1, x_2 ) ! Out
! Description:
! For the Luhar closure, this subroutine takes Skx, xm, xp2, and mixt_frac,
! big_m, and small_m (calculated in calc_Luhar_params) as input and outputs
! the PDF component means and variances of a variable x in the joint-PDF
! according to Luhar et al. (1996). This code was written using the
! equations and nomenclature of Larson et al. (2002) Appendix section e.
! References:
! Vincent E. Larson, Jean-Christophe Golaz, and William R. Cotton, 2002:
! Small-Scale and Mesoscale Variability in Cloudy Boundary Layers: Joint
! Probability Density Functions. J. Atmos. Sci., 59, 3519–3539.
!-----------------------------------------------------------------------
use constants_clubb, only: &
one, & ! Constant(s)
zero
use clubb_precision, only: &
core_rknd ! Precision
implicit none
! Input Variables
integer, intent(in) :: &
nz ! Number of vertical levels
real( kind = core_rknd ), dimension(nz), intent(in) :: &
xm, & ! Mean (overall) of x, <x> [(x units)]
xp2, & ! Variance (overall) of x, <x'^2> [(x units)^2]
mixt_frac, & ! Mixture fraction [-]
small_m, & ! Luhar's small m [-]
wpxp ! Covariance of w and x, <w'x'> [m/s (x units)]
real( kind = core_rknd ), intent(in) :: &
x_tol_sqd ! Tolerance value of x squared [(x units)^2]
! Output Variables
real( kind = core_rknd ), dimension(nz), intent(out) :: &
sigma_sqd_x_1, & ! Normalized width parameter of x (1st PDF component) [-]
sigma_sqd_x_2, & ! Normalized width parameter of x (1st PDF component) [-]
varnce_x_1, & ! Variance of x (1st PDF component) [(x units)^2]
varnce_x_2, & ! Variance of x (2nd PDF component) [(x units)^2]
x_1_n, & ! Normalized mean of x (1st PDF component) [-]
x_2_n, & ! Normalized mean of x (2nd PDF component) [-]
x_1, & ! Mean of x (1st PDF component) [(x units)]
x_2 ! Mean of x (2nd PDF component) [(x units)]
! Local Variables
real( kind = core_rknd), dimension(nz) :: &
sqrt_xp2, & ! Square root of the variance of x [(x units)]
sgn_wpxp ! Sgn( <w'x'> ); 1 when <w'x'> >= 0 or -1 when <w'x'> < 0 [-]
where ( xp2 > x_tol_sqd )
! Width (standard deviation) parameters are non-zero
! Calculate sgn( <w'x'> ).
where ( wpxp >= zero )
sgn_wpxp = one
elsewhere ! <w'x'> < 0
sgn_wpxp = -one
endwhere ! <w'x'> >= 0
! Calculate the square root of the overall variance of x.
sqrt_xp2 = sqrt( xp2 )
! Normalized width parameter of x in the 1st PDF component.
sigma_sqd_x_1 &
= ( one - mixt_frac ) / ( mixt_frac * ( one + small_m**2 ) )
! The variance of x in the 1st PDF component.
varnce_x_1 = sigma_sqd_x_1 * xp2
! Normalized width parameter of x in the 2nd PDF component.
sigma_sqd_x_2 &
= mixt_frac / ( ( one - mixt_frac ) * ( one + small_m**2 ) )
! The variance of x in the 2nd PDF component.
varnce_x_2 = sigma_sqd_x_2 * xp2
! Normalized mean of x in the 1st PDF component.
x_1_n = sgn_wpxp * small_m * sqrt( sigma_sqd_x_1 )
! Normalized mean of x in the 2nd PDF component.
x_2_n = -sgn_wpxp * small_m * sqrt( sigma_sqd_x_2 )
! The mean of x in the 1st PDF component.
x_1 = xm + sqrt_xp2 * x_1_n
! The mean of x in the 2nd PDF component.
x_2 = xm + sqrt_xp2 * x_2_n
elsewhere
! Variable x doesn't vary.
sigma_sqd_x_1 = ( one / ( one + small_m**2 ) )
sigma_sqd_x_2 = ( one / ( one + small_m**2 ) )
varnce_x_1 = zero
varnce_x_2 = zero
x_1_n = sgn_wpxp * small_m * sqrt( sigma_sqd_x_1 )
x_2_n = -sgn_wpxp * small_m * sqrt( sigma_sqd_x_2 )
x_1 = xm
x_2 = xm
endwhere ! Widths non-zero
return
end subroutine close_Luhar_pdf
!=============================================================================
subroutine ADG1_ADG2_responder_params( xm, xp2, wp2, sqrt_wp2, & ! In
wpxp, w_1_n, w_2_n, mixt_frac, & ! In
sigma_sqd_w, x_tol, & ! In
x_1, x_2, varnce_x_1, & ! Out
varnce_x_2, alpha_x ) ! Out
! Description:
! Calculates the PDF component means and PDF component variances for thl,
! rt, or sclr when either the ADG1 PDF or ADG2 PDF are used.
!
! The normalized variance for thl, rt, and sclr for "plume" 1 is:
!
! { 1 - [1/(1-sigma_sqd_w)]*[ <w'x'>^2 / (<w'^2> * <x'^2>) ] / mixt_frac }
! * { (1/3)*beta + mixt_frac*( 1 - (2/3)*beta ) };
!
! where "x" stands for thl, rt, or sclr; "mixt_frac" is the weight of
! Gaussian "plume" 1, and 0 <= beta <= 3.
!
! The factor { (1/3)*beta + mixt_frac*( 1 - (2/3)*beta ) } does not depend
! on which varable "x" stands for. The factor is multiplied by 2 and
! defined as width_factor_1.
!
! The factor
! { 1 - [1/(1-sigma_sqd_w)]*[ (w'x')^2 / (w'^2 * x'^2) ] / mixt_frac }
! depends on which variable "x" stands for. It is multiplied by 1/2 and
! defined as alpha_x, where "x" stands for thl, rt, or sclr.
! Vince Larson added a dimensionless factor so that the
! width of plumes in theta_l, rt can vary.
! beta is a constant defined in module parameters_tunable
! Set 0<beta<3.
! beta=1.5_core_rknd recovers Chris Golaz' simplified formula.
! 3 Nov 2003
! References:
use grid_class, only: &
gr ! Variable type(s)
use constants_clubb, only: &
two, & ! Constant(s)
one, &
two_thirds, &
one_half, &
zero, &
zero_threshold, &
w_tol_sqd
use parameters_tunable, only: &
beta ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Input Variables
real ( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
xm, & ! Mean of x (overall) [units vary]