-
Notifications
You must be signed in to change notification settings - Fork 12
/
parameters_tunable.F90
2099 lines (1788 loc) · 81.1 KB
/
parameters_tunable.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 parameters_tunable
! Description:
! This module contains tunable model parameters. The purpose of the module is to make it
! easier for the clubb_tuner code to use the params vector without "knowing" any information
! about the individual parameters contained in the vector itself. It makes it easier to add
! new parameters to be tuned for, but does not make the CLUBB_core code itself any simpler.
! The parameters within the vector do not need to be the same variables used in the rest of
! CLUBB_core (see for e.g. nu1_vert_res_dep or lmin_coef).
! The parameters in the params vector only need to be those parameters for which we're not
! sure the correct value and we'd like to tune for.
!
! References:
! None
!
! Notes:
! To make it easier to verify of code correctness, please keep the omp threadprivate
! directives just after the variable declaration. All parameters in this
! module should be declared threadprivate because of the CLUBB tuner.
!-----------------------------------------------------------------------
use parameter_indices, only: nparams ! Variable(s)
use grid_class, only: gr ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Default to private
private
public :: setup_parameters, read_parameters, read_param_max, &
get_parameters, adj_low_res_nu, cleanup_nu
! NOTE: In CLUBB standalone, as well as some host models, the hardcoded
! default values of some or all of the parameters below have no effect,
! as the values are simply read in using a namelist or set in host model
! specific code.
! Model constant parameters
real( kind = core_rknd ), public :: &
C1 = 1.000000_core_rknd, & ! Low Skewness in C1 Skw. Function [-]
C1b = 1.000000_core_rknd, & ! High Skewness in C1 Skw. Function [-]
C1c = 1.000000_core_rknd, & ! Degree of Slope of C1 Skw. Function [-]
C2 = 1.300000_core_rknd, & ! Low Skewness in C2 Skw. Function [-]
C2rt = 1.000000_core_rknd, & ! C2 coef. for the rtp2_dp1 term [-]
C2thl = 1.000000_core_rknd, & ! C2 coef. for the thlp2_dp1 term [-]
C2rtthl = 2.000000_core_rknd, & ! C2 coef. for the rtpthlp_dp1 term [-]
C2b = 1.300000_core_rknd, & ! High Skewness in C2 Skw. Function [-]
C2c = 5.000000_core_rknd, & ! Degree of Slope of C2 Skw. Function [-]
C4 = 5.200000_core_rknd, & ! Used only when l_tke_aniso is true [-]
C5 = 0.300000_core_rknd, & ! Coef. in pressure terms: w'^2 eqn [-]
C6rt = 4.000000_core_rknd, & ! Low Skewness in C6rt Skw. Function [-]
C6rtb = 6.000000_core_rknd, & ! High Skewness in C6rt Skw. Function [-]
C6rtc = 1.000000_core_rknd, & ! Degree of Slope of C6rt Skw. Fnct. [-]
C6thl = 4.000000_core_rknd, & ! Low Skewness in C6thl Skw. Function [-]
C6thlb = 6.000000_core_rknd, & ! High Skewness in C6thl Skw. Fnct. [-]
C6thlc = 1.000000_core_rknd, & ! Degree of Slope of C6thl Skw. Fnct. [-]
C7 = 0.500000_core_rknd, & ! Low Skewness in C7 Skw. Function [-]
C7b = 0.800000_core_rknd, & ! High Skewness in C7 Skw. Function [-]
C7c = 0.500000_core_rknd, & ! Degree of Slope of C7 Skw. Function [-]
C8 = 3.000000_core_rknd, & ! Coef. #1 in C8 Skewness Equation [-]
C8b = 0.020000_core_rknd, & ! Coef. #2 in C8 Skewness Equation [-]
C10 = 3.300000_core_rknd, & ! Currently Not Used in the Model [-]
C11 = 0.500000_core_rknd, & ! Low Skewness in C11 Skw. Function [-]
C11b = 0.500000_core_rknd, & ! High Skewness in C11 Skw. Function [-]
C11c = 0.500000_core_rknd, & ! Degree of Slope of C11 Skw. Fnct. [-]
C12 = 1.000000_core_rknd, & ! Constant in w'^3 Crank-Nich. diff. [-]
C13 = 0.100000_core_rknd, & ! Not currently used in model [-]
C14 = 1.000000_core_rknd, & ! Constant for u'^2 and v'^2 terms [-]
C15 = 0.4_core_rknd, & ! Coefficient for the wp3_bp2 term [-]
C_wp2_splat = 0.0_core_rknd ! Coefficient for gustiness near ground [-]
!$omp threadprivate(C1, C1b, C1c, C2, C2b, C2c, &
!$omp C2rt, C2thl, C2rtthl, C4, C5, C6rt, C6rtb, C6rtc, &
!$omp C6thl, C6thlb, C6thlc, &
!$omp C7, C7b, C7c, C8, C8b, C10, C11, C11b, C11c, C12, &
!$omp C13, C14, C15, C_wp2_splat)
real( kind = core_rknd ), public :: &
C6rt_Lscale0 = 14.0_core_rknd, & ! Damp C6rt as a fnct. of Lscale [-]
C6thl_Lscale0 = 14.0_core_rknd, & ! Damp C6thl as a fnct. of Lscale [-]
C7_Lscale0 = 0.8500000_core_rknd, & ! Damp C7 as a fnct. of Lscale [-]
wpxp_L_thresh = 60.0_core_rknd ! Lscale threshold: damp C6 & C7 [m]
!$omp threadprivate(C6rt_Lscale0, C6thl_Lscale0, C7_Lscale0, wpxp_L_thresh)
! Note: DD 1987 is Duynkerke & Driedonks (1987).
real( kind = core_rknd ), public :: &
c_K = 0.200000_core_rknd, & ! Constant C_mu^(1/4) in DD 1987 [m^2/s]
c_K1 = 0.750000_core_rknd, & ! Coef. of Eddy Diffusion: wp2 [m^2/s]
c_K2 = 0.125000_core_rknd, & ! Coef. of Eddy Diffusion: xp2 [m^2/s]
c_K6 = 0.375000_core_rknd, & ! Coef. of Eddy Diffusion: wpxp [m^2/s]
c_K8 = 1.250000_core_rknd, & ! Coef. of Eddy Diffusion: wp3 [m^2/s]
c_K9 = 0.250000_core_rknd, & ! Coef. of Eddy Diff.: up2/vp2 [m^2/s]
c_K_hm = 0.750000_core_rknd, & ! Coef. of Eddy Diffusion: hmm [m^2/s]
c_K_hmb = 0.10000_core_rknd, & ! Coef. of Non-Local Factor, Eddy Diffusion: hmm [m^2/s]
K_hm_min_coef = 0.10000_core_rknd,& ! Min. of Non-Local Factor, Eddy Diffusion: hmm [m^2/s]
gamma_coef = 0.320000_core_rknd, & ! Low Skw.: gamma coef. Skw. Fnct. [-]
gamma_coefb = 0.320000_core_rknd, & ! High Skw.: gamma coef. Skw. Fnct. [-]
gamma_coefc = 5.000000_core_rknd, & ! Deg. Slope: gamma coef. Skw. Fnct. [-]
mu = 1.000E-3_core_rknd, & ! Fract entrain rate per unit alt [1/m]
mult_coef = 0.500000_core_rknd, & ! Coef. applied to log(avg dz/thresh)[-]
taumin = 90.00000_core_rknd, & ! Min. allow. value: time-scale tau [s]
taumax = 3600.000_core_rknd, & ! Max. allow. value: time-scale tau [s]
lmin = 20.00000_core_rknd ! Min. value for the length scale [m]
!$omp threadprivate(c_K, c_K1, c_K2, c_K6, &
!$omp c_K8, c_K9, c_K_hm, c_K_hmb, K_hm_min_coef, gamma_coef, gamma_coefb, gamma_coefc, &
!$omp mu, mult_coef, taumin, taumax, lmin)
real( kind = core_rknd ), public :: &
Lscale_mu_coef = 2.0_core_rknd, & ! Coef perturb mu: av calc Lscale [-]
Lscale_pert_coef = 0.1_core_rknd ! Coef pert thlm/rtm: av calc Lscale [-]
!$omp threadprivate(Lscale_mu_coef, Lscale_pert_coef)
real( kind = core_rknd ), public :: &
alpha_corr = 0.15_core_rknd ! Coef. for the corr. diagnosis algorithm [-]
!$omp threadprivate(alpha_corr)
real( kind = core_rknd ), private :: &
nu1 = 20.00000_core_rknd, & ! Bg. Coef. Eddy Diffusion: wp2 [m^2/s]
nu2 = 5.000000_core_rknd, & ! Bg. Coef. Eddy Diffusion: xp2 [m^2/s]
nu6 = 5.000000_core_rknd, & ! Bg. Coef. Eddy Diffusion: wpxp [m^2/s]
nu8 = 20.00000_core_rknd, & ! Bg. Coef. Eddy Diffusion: wp3 [m^2/s]
nu9 = 20.00000_core_rknd, & ! Bg. Coef. Eddy Diffusion: up2/vp2 [m^2/s]
nu10 = 0.000000_core_rknd, & ! Bg. Coef. Eddy Diffusion: edsclrm [m^2/s]
nu_hm = 1.500000_core_rknd ! Bg. Coef. Eddy Diffusion: hmm [m^2/s]
!$omp threadprivate(nu1, nu2, nu6, nu8, nu9, nu10, nu_hm)
real( kind = core_rknd ), public, allocatable, dimension(:) :: &
nu1_vert_res_dep, & ! Background Coef. of Eddy Diffusion: wp2 [m^2/s]
nu2_vert_res_dep, & ! Background Coef. of Eddy Diffusion: xp2 [m^2/s]
nu6_vert_res_dep, & ! Background Coef. of Eddy Diffusion: wpxp [m^2/s]
nu8_vert_res_dep, & ! Background Coef. of Eddy Diffusion: wp3 [m^2/s]
nu9_vert_res_dep, & ! Background Coef. of Eddy Diffusion: up2/vp2 [m^2/s]
nu10_vert_res_dep, & ! Background Coef. of Eddy Diffusion: edsclrm [m^2/s]
nu_hm_vert_res_dep ! Background Coef. of Eddy Diffusion: hydromet [m^2/s]
!$omp threadprivate(nu1_vert_res_dep, nu2_vert_res_dep, nu6_vert_res_dep, &
!$omp nu8_vert_res_dep, nu9_vert_res_dep, nu10_vert_res_dep, nu_hm_vert_res_dep)
! Vince Larson added a constant to set plume widths for theta_l and rt
! beta should vary between 0 and 3.
real( kind = core_rknd ), public :: &
beta = 2.400000_core_rknd ! Beta coefficient [-]
!$omp threadprivate(beta)
real( kind = core_rknd ), private :: &
lmin_coef = 0.500000_core_rknd ! Coefficient of lmin [-]
!$omp threadprivate(lmin_coef)
real( kind = core_rknd ), public :: &
Skw_max_mag = 10.0_core_rknd ! Max magnitude of skewness [-]
!$omp threadprivate(Skw_max_mag)
real( kind = core_rknd ), public :: &
C_invrs_tau_bkgnd = 1.0_core_rknd, & !
C_invrs_tau_sfc = 0.1_core_rknd, & !
C_invrs_tau_shear = 0.02_core_rknd, & !
C_invrs_tau_N2 = 0.1_core_rknd, & !
C_invrs_tau_N2_wp2 = 0.2_core_rknd, & !
C_invrs_tau_N2_xp2 = 0.2_core_rknd, & !
C_invrs_tau_N2_wpxp = 0.0_core_rknd, & !
C_invrs_tau_N2_clear_wp3 = 0.0_core_rknd
!$omp threadprivate(C_invrs_tau_bkgnd,C_invrs_tau_sfc)
!$omp threadprivate(C_invrs_tau_shear,C_invrs_tau_N2)
!$omp threadprivate(C_invrs_tau_N2_wp2,C_invrs_tau_N2_xp2)
!$omp threadprivate(C_invrs_tau_N2_wpxp)
!$omp threadprivate(C_invrs_tau_N2_clear_wp3)
! Parameters for the new PDF (w, rt, and theta-l).
!
! Brian Griffin added a tunable parameter for the PDF of w,
! slope_coef_spread_DG_means_w, to increase or decrease the spread between the
! two PDF component means of w. When the value of this slope parameter is
! larger, F_w is smaller and the PDF component means of w are closer together.
! Valid values are slope_coef_spread_DG_means_w > 0.
!
! A second parameter for the PDF of w, pdf_component_stdev_factor_w, is used
! to adjust the standard deviations of the 1st PDF component against the 2nd
! PDF component for w. This parameter is related to zeta_w, where:
!
! 1 + zeta_w = ( mixt_frac * sigma_w_1^2 )
! / ( ( 1 - mixt_frac ) * sigma_w_2^2 );
!
! The pdf_component_stdev_factor_w is set such that:
!
! pdf_component_stdev_factor_w = zeta_w + 1.
!
! Valid values are pdf_component_stdev_factor_w > 0.
!
! The parameter for the PDF of rt is coef_spread_DG_means_rt. Valid values
! are 0 <= coef_spread_DG_means_rt < 1. When coef_spread_DG_means_rt
! approaches 0, F_rt approaches min_F_rt, and the two PDF component means
! become closer together. When coef_spread_DG_means_rt approaches 1, F_rt
! approaches max_F_rt, and the two PDF component means are spread farther
! apart.
!
! The parameter for the PDF of theta-l is coef_spread_DG_means_thl.
! Valid values are 0 <= coef_spread_DG_means_thl < 1. When
! coef_spread_DG_means_thl approaches 0, F_thl approaches min_F_thl, and the
! two PDF component means become closer together. When
! coef_spread_DG_means_thl approaches 1, F_thl approaches max_F_thl, and the
! two PDF component means are spread farther apart.
real( kind = core_rknd ), public :: &
! Slope coefficient for the spread between the PDF component means of w.
slope_coef_spread_DG_means_w = 21.0_core_rknd, &
! Parameter to adjust the PDF component standard deviations of w.
pdf_component_stdev_factor_w = 6.5_core_rknd, &
! Coefficient for the spread between the PDF component means of rt.
coef_spread_DG_means_rt = 0.8_core_rknd, &
! Coefficient for the spread between the PDF component means of thl.
coef_spread_DG_means_thl = 0.8_core_rknd
!$omp threadprivate( slope_coef_spread_DG_means_w, &
!$omp pdf_component_stdev_factor_w, &
!$omp coef_spread_DG_means_rt, &
!$omp coef_spread_DG_means_thl )
! Parameters for the hydrometeor portion of the PDF.
!
! Brian Griffin added a parameter for hydrometeors, omicron, to increase the
! standard deviation of each component and decrease the spread between the
! component means as the value of omicron inreases. Valid value are
! 0 < omicron <= 1.
! A second parameter for hydrometeors, zeta, increases the standard deviation
! of component 1 at the expense of the standard deviation of component 2 when
! the value of zeta > 0 (and increasingly so as zeta increases). Valid values
! are zeta > -1.
real( kind = core_rknd ), public :: &
omicron = 0.8_core_rknd, & ! Hydromet width/spread-of-means param [-]
zeta_vrnce_rat = 0.0_core_rknd ! Ratio sigma^2/mu^2 comp. 1 / comp. 2 [-]
!$omp threadprivate( omicron, zeta_vrnce_rat )
real( kind = core_rknd ), public :: &
! ratio mixt_frac*precip_frac_1/precip_frac (precip_frac_calc_type=2) [-]
upsilon_precip_frac_rat = 0.9_core_rknd, &
! Intensity of stability correction applied to C1 and C6 [-]
lambda0_stability_coef = 0.03_core_rknd
!$omp threadprivate( upsilon_precip_frac_rat, lambda0_stability_coef)
! Factor to decrease sensitivity in the denominator of Skw calculation
real( kind = core_rknd ), public :: &
Skw_denom_coef = 4.0_core_rknd
!$omp threadprivate( Skw_denom_coef )
! Momentum coefficient of Kh_zm
real( kind = core_rknd ), public :: &
c_K10 = 1.0_core_rknd
!$omp threadprivate( c_K10 )
! Thermodynamic coefficient of Kh_zm
real( kind = core_rknd ), public :: &
c_K10h = 1.0_core_rknd
!$omp threadprivate( c_K10h )
real( kind = core_rknd ), public :: &
thlp2_rad_coef = 1.0_core_rknd, & ! Coefficient of thlp2_rad [-]
thlp2_rad_cloud_frac_thresh = 0.1_core_rknd ! Minimum cloud fraction for computation
! of thlp2_rad [-]
!$omp threadprivate( thlp2_rad_coef, thlp2_rad_cloud_frac_thresh )
real( kind = core_rknd ), public :: &
up2_vp2_factor = 2.0_core_rknd ! Coefficients of up2 and vp2 [-]
!$omp threadprivate( up2_vp2_factor )
! Since we lack a devious way to do this just once, this namelist
! must be changed as well when a new parameter is added.
namelist /clubb_params_nl/ &
C1, C1b, C1c, C2, C2b, C2c, &
C2rt, C2thl, C2rtthl, C4, C5, &
C6rt, C6rtb, C6rtc, C6thl, C6thlb, C6thlc, &
C7, C7b, C7c, C8, C8b, C10, C11, C11b, C11c, &
C12, C13, C14, C15, C_wp2_splat, &
C6rt_Lscale0, C6thl_Lscale0, &
C7_Lscale0, wpxp_L_thresh, c_K, c_K1, nu1, c_K2, nu2, &
c_K6, nu6, c_K8, nu8, c_K9, nu9, nu10, &
c_K_hm, c_K_hmb, K_hm_min_coef, nu_hm, &
slope_coef_spread_DG_means_w, pdf_component_stdev_factor_w, &
coef_spread_DG_means_rt, coef_spread_DG_means_thl, &
beta, gamma_coef, gamma_coefb, gamma_coefc, lmin_coef, &
omicron, zeta_vrnce_rat, upsilon_precip_frac_rat, &
lambda0_stability_coef, mult_coef, taumin, taumax, mu, Lscale_mu_coef, &
Lscale_pert_coef, alpha_corr, Skw_denom_coef, c_K10, c_K10h, &
thlp2_rad_coef, thlp2_rad_cloud_frac_thresh, up2_vp2_factor, Skw_max_mag, &
C_invrs_tau_bkgnd, C_invrs_tau_sfc, C_invrs_tau_shear, C_invrs_tau_N2, &
C_invrs_tau_N2_wp2, C_invrs_tau_N2_xp2, &
C_invrs_tau_N2_wpxp, C_invrs_tau_N2_clear_wp3
! These are referenced together often enough that it made sense to
! make a list of them. Note that lmin_coef is the input parameter,
! while the actual lmin model constant is computed from this.
!***************************************************************
! ***** IMPORTANT *****
! If you change the order of the parameters in the parameter_indices,
! you will need to change the order of this list as well or the
! tuner will break!
! ***** IMPORTANT *****
!***************************************************************
character(len=28), dimension(nparams), parameter, public :: &
params_list = &
(/"C1 ", "C1b ", &
"C1c ", "C2 ", &
"C2b ", "C2c ", &
"C2rt ", "C2thl ", &
"C2rtthl ", "C4 ", &
"C5 ", "C6rt ", &
"C6rtb ", "C6rtc ", &
"C6thl ", "C6thlb ", &
"C6thlc ", "C7 ", &
"C7b ", "C7c ", &
"C8 ", "C8b ", &
"C10 ", "C11 ", &
"C11b ", "C11c ", &
"C12 ", "C13 ", &
"C14 ", "C15 ", &
"C_wp2_splat ", &
"C6rt_Lscale0 ", "C6thl_Lscale0 ", &
"C7_Lscale0 ", "wpxp_L_thresh ", &
"c_K ", "c_K1 ", &
"nu1 ", "c_K2 ", &
"nu2 ", "c_K6 ", &
"nu6 ", "c_K8 ", &
"nu8 ", "c_K9 ", &
"nu9 ", "nu10 ", &
"c_K_hm ", "c_K_hmb ", &
"K_hm_min_coef ", "nu_hm ", &
"slope_coef_spread_DG_means_w", "pdf_component_stdev_factor_w", &
"coef_spread_DG_means_rt ", "coef_spread_DG_means_thl ", &
"gamma_coef ", "gamma_coefb ", &
"gamma_coefc ", "mu ", &
"beta ", "lmin_coef ", &
"omicron ", "zeta_vrnce_rat ", &
"upsilon_precip_frac_rat ", "lambda0_stability_coef ", &
"mult_coef ", "taumin ", &
"taumax ", "Lscale_mu_coef ", &
"Lscale_pert_coef ", "alpha_corr ", &
"Skw_denom_coef ", "c_K10 ", &
"c_K10h ", "thlp2_rad_coef ", &
"thlp2_rad_cloud_frac_thresh ", "up2_vp2_factor ", &
"Skw_max_mag ", "C_invrs_tau_bkgnd ", &
"C_invrs_tau_sfc ", "C_invrs_tau_shear ", &
"C_invrs_tau_N2 ", "C_invrs_tau_N2_wp2 ", &
"C_invrs_tau_N2_xp2 ", "C_invrs_tau_N2_wpxp ", &
"C_invrs_tau_N2_clear_wp3 " /)
real( kind = core_rknd ), parameter, private :: &
init_value = -999._core_rknd ! Initial value for the parameters, used to detect missing values
#ifdef E3SM
public :: clubb_param_readnl
! The parameters below have the same meaning as those without prefix 'clubb_'
! They can be specified via namelist to ease parameter tuning
! If adding more parameters for tuning via namelist, need to insert blocks
! accordingly in subroutines read_parameters and clubb_paramm_readnl
! (including updating list of nml variables)
real( kind = core_rknd ) :: &
clubb_C1, &
clubb_C1b, &
clubb_C1c, &
clubb_C2rt, &
clubb_C2thl, &
clubb_C2rtthl, &
clubb_C4, &
clubb_C5, &
clubb_C6rt, &
clubb_C6rtb, &
clubb_C6rtc, &
clubb_C6thlb, &
clubb_C6thlc, &
clubb_C7, &
clubb_C7b, &
clubb_C8, &
clubb_C11, &
clubb_C11b, &
clubb_C11c, &
clubb_C14, &
clubb_C15, &
clubb_beta, &
clubb_gamma_coef, &
clubb_gamma_coefb, &
clubb_gamma_coefc, &
clubb_mu, &
clubb_nu1, &
clubb_nu2, &
clubb_c_K2, &
clubb_c_K10, &
clubb_wpxp_L_thresh, &
clubb_C_invrs_tau_bkgnd, &
clubb_C_invrs_tau_sfc, &
clubb_C_invrs_tau_shear, &
clubb_C_invrs_tau_N2, &
clubb_C_invrs_tau_N2_wp2, &
clubb_C_invrs_tau_N2_xp2, &
clubb_C_invrs_tau_N2_wpxp, &
clubb_C_invrs_tau_N2_clear_wp3,&
clubb_C_wp2_splat
#endif /*E3SM*/
contains
!=============================================================================
subroutine setup_parameters &
( deltaz, params, nzmax, &
grid_type, momentum_heights, thermodynamic_heights, &
l_prescribed_avg_deltaz, &
err_code_out )
! Description:
! Subroutine to setup model parameters
! References:
! None
!-----------------------------------------------------------------------
use constants_clubb, only: &
three, & ! Variable(s)
one, &
zero, &
fstderr
use model_flags, only: &
l_clip_semi_implicit ! Variable(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
use error_code, only: &
err_code, & ! Error Indicator
clubb_fatal_error ! Constant
use parameter_indices, only: &
izeta_vrnce_rat
implicit none
! Constant Parameters
real( kind = core_rknd ), parameter :: &
lmin_deltaz = 40.0_core_rknd ! Fixed value for minimum value for the length scale.
! Input Variables
real( kind = core_rknd ), intent(in) :: &
deltaz ! Change per height level [m]
real( kind = core_rknd ), intent(in), dimension(nparams) :: &
params ! Tuneable model parameters [-]
! Grid definition
integer, intent(in) :: nzmax ! Vertical grid levels [#]
! If CLUBB is running on its own, this option determines
! if it is using:
! 1) an evenly-spaced grid,
! 2) a stretched (unevenly-spaced) grid entered on the
! thermodynamic grid levels (with momentum levels set
! halfway between thermodynamic levels), or
! 3) a stretched (unevenly-spaced) grid entered on the
! momentum grid levels (with thermodynamic levels set
! halfway between momentum levels).
integer, intent(in) :: grid_type
! If the CLUBB parameterization is implemented in a host model,
! it needs to use the host model's momentum level altitudes
! and thermodynamic level altitudes.
! If the CLUBB model is running by itself, but is using a
! stretched grid entered on thermodynamic levels (grid_type = 2),
! it needs to use the thermodynamic level altitudes as input.
! If the CLUBB model is running by itself, but is using a
! stretched grid entered on momentum levels (grid_type = 3),
! it needs to use the momentum level altitudes as input.
real( kind = core_rknd ), intent(in), dimension(nzmax) :: &
momentum_heights, & ! Momentum level altitudes (input) [m]
thermodynamic_heights ! Thermodynamic level altitudes (input) [m]
logical, intent(in) :: &
l_prescribed_avg_deltaz ! used in adj_low_res_nu. If .true., avg_deltaz = deltaz
integer, intent(out) :: &
err_code_out ! Error code indicator
integer :: k ! loop variable
!-------------------- Begin code --------------------
! Ensure all variables are greater than 0, and zeta_vrnce_rat is greater than -1
do k = 1, nparams
if ( k /= izeta_vrnce_rat .and. params(k) < zero ) then
write(fstderr,*) params_list(k), " = ", params(k)
write(fstderr,*) params_list(k), " must satisfy 0.0 <= ", params_list(k)
err_code = clubb_fatal_error
else if ( params(k) < -one ) then
write(fstderr,*) "zeta_vrnce_rat = ", zeta_vrnce_rat
write(fstderr,*) "zeta_vrnce_rat must satisfy -1.0 <= zeta_vrnce_rat"
err_code = clubb_fatal_error
end if
end do
call unpack_parameters &
( params, &
C1, C1b, C1c, C2, C2b, C2c, C2rt, C2thl, C2rtthl, &
C4, C5, C6rt, C6rtb, C6rtc, C6thl, C6thlb, C6thlc, &
C7, C7b, C7c, C8, C8b, C10, &
C11, C11b, C11c, C12, C13, C14, C15, C_wp2_splat, &
C6rt_Lscale0, C6thl_Lscale0, C7_Lscale0, wpxp_L_thresh, &
c_K, c_K1, nu1, c_K2, nu2, c_K6, nu6, c_K8, nu8, &
c_K9, nu9, nu10, c_K_hm, c_K_hmb, K_hm_min_coef, nu_hm, &
slope_coef_spread_DG_means_w, pdf_component_stdev_factor_w, &
coef_spread_DG_means_rt, coef_spread_DG_means_thl, &
gamma_coef, gamma_coefb, gamma_coefc, mu, beta, lmin_coef, &
omicron, zeta_vrnce_rat, upsilon_precip_frac_rat, &
lambda0_stability_coef, mult_coef, taumin, taumax, &
Lscale_mu_coef, Lscale_pert_coef, alpha_corr, &
Skw_denom_coef, c_K10, c_K10h, thlp2_rad_coef, &
thlp2_rad_cloud_frac_thresh, up2_vp2_factor, Skw_max_mag, &
C_invrs_tau_bkgnd, C_invrs_tau_sfc, &
C_invrs_tau_shear, C_invrs_tau_N2, &
C_invrs_tau_N2_wp2, C_invrs_tau_N2_xp2, &
C_invrs_tau_N2_wpxp, C_invrs_tau_N2_clear_wp3 )
! It was decided after some experimentation, that the best
! way to produce grid independent results is to set lmin to be
! some fixed value. -dschanen 21 May 2007
!lmin = lmin_coef * deltaz ! Old
lmin = lmin_coef * lmin_deltaz ! New fixed value
! ### Adjust Constant Diffusivity Coefficients Based On Grid Spacing ###
call adj_low_res_nu &
( nzmax, grid_type, deltaz, & ! Intent(in)
momentum_heights, thermodynamic_heights, & ! Intent(in)
l_prescribed_avg_deltaz ) ! Intent(in)
if ( beta < zero .or. beta > three ) then
! Constraints on beta
write(fstderr,*) "beta = ", beta
write(fstderr,*) "beta cannot be < 0 or > 3"
err_code = clubb_fatal_error
endif ! beta < 0 or beta > 3
if ( slope_coef_spread_DG_means_w <= zero ) then
! Constraint on slope_coef_spread_DG_means_w
write(fstderr,*) "slope_coef_spread_DG_means_w = ", &
slope_coef_spread_DG_means_w
write(fstderr,*) "slope_coef_spread_DG_means_w cannot be <= 0"
err_code = clubb_fatal_error
endif ! slope_coef_spread_DG_means_w <= 0
if ( pdf_component_stdev_factor_w <= zero ) then
! Constraint on pdf_component_stdev_factor_w
write(fstderr,*) "pdf_component_stdev_factor_w = ", &
pdf_component_stdev_factor_w
write(fstderr,*) "pdf_component_stdev_factor_w cannot be <= 0"
err_code = clubb_fatal_error
endif ! pdf_component_stdev_factor_w <= 0
if ( coef_spread_DG_means_rt < zero &
.or. coef_spread_DG_means_rt >= one ) then
! Constraint on coef_spread_DG_means_rt
write(fstderr,*) "coef_spread_DG_means_rt = ", coef_spread_DG_means_rt
write(fstderr,*) "coef_spread_DG_means_rt cannot be < 0 or >= 1"
err_code = clubb_fatal_error
endif ! coef_spread_DG_means_rt < 0 or coef_spread_DG_means_rt >= 1
if ( coef_spread_DG_means_thl < zero &
.or. coef_spread_DG_means_thl >= one ) then
! Constraint on coef_spread_DG_means_thl
write(fstderr,*) "coef_spread_DG_means_thl = ", coef_spread_DG_means_thl
write(fstderr,*) "coef_spread_DG_means_thl cannot be < 0 or >= 1"
err_code = clubb_fatal_error
endif ! coef_spread_DG_means_thl < 0 or coef_spread_DG_means_thl >= 1
if ( omicron <= zero .or. omicron > one ) then
! Constraints on omicron
write(fstderr,*) "omicron = ", omicron
write(fstderr,*) "omicron cannot be <= 0 or > 1"
err_code = clubb_fatal_error
endif ! omicron <= 0 or omicron > 1
if ( zeta_vrnce_rat <= -one ) then
! Constraints on zeta_vrnce_rat
write(fstderr,*) "zeta_vrnce_rat = ", zeta_vrnce_rat
write(fstderr,*) "zeta_vrnce_rat cannot be <= -1"
err_code = clubb_fatal_error
endif ! zeta_vrnce_rat <= -1
if ( upsilon_precip_frac_rat < zero &
.or. upsilon_precip_frac_rat > one ) then
! Constraints on upsilon_precip_frac_rat
write(fstderr,*) "upsilon_precip_frac_rat = ", upsilon_precip_frac_rat
write(fstderr,*) "upsilon_precip_frac_rat cannot be < 0 or > 1"
err_code = clubb_fatal_error
endif ! upsilon_precip_frac_rat < 0 or upsilon_precip_frac_rat > 1
if ( mu < zero ) then
! Constraints on entrainment rate, mu.
write(fstderr,*) "mu = ", mu
write(fstderr,*) "mu cannot be < 0"
err_code = clubb_fatal_error
endif ! mu < 0.0
if ( lmin < 4.0_core_rknd ) then
! Constraints on mixing length
write(fstderr,*) "lmin = ", lmin
write(fstderr,*) "lmin is < 4.0_core_rknd"
err_code = clubb_fatal_error
endif ! lmin < 4.0
if ( .not. l_clip_semi_implicit ) then
! When l_clip_semi_implicit is set to false (the current default),
! the C6rt parameters must be set equal to the C6thl parameters.
! Otherwise, the wpthlp pr1 term will be calculated inconsistently.
if ( C6rt /= C6thl ) then
write(fstderr,*) "C6rt = ", C6rt
write(fstderr,*) "C6thl = ", C6thl
write(fstderr,*) "C6rt and C6thl must be equal when" &
// " l_clip_semi_implicit is turned off (default)."
err_code = clubb_fatal_error
endif ! C6rt /= C6thl
if ( C6rtb /= C6thlb ) then
write(fstderr,*) "C6rtb = ", C6rtb
write(fstderr,*) "C6thlb = ", C6thlb
write(fstderr,*) "C6rtb and C6thlb must be equal when" &
// " l_clip_semi_implicit is turned off (default)."
err_code = clubb_fatal_error
endif ! C6rtb /= C6thlb
if ( C6rtc /= C6thlc ) then
write(fstderr,*) "C6rtc = ", C6rtc
write(fstderr,*) "C6thlc = ", C6thlc
write(fstderr,*) "C6rtc and C6thlc must be equal when" &
// " l_clip_semi_implicit is turned off (default)."
err_code = clubb_fatal_error
endif ! C6rtc /= C6thlc
if ( C6rt_Lscale0 /= C6thl_Lscale0 ) then
write(fstderr,*) "C6rt_Lscale0 = ", C6rt_Lscale0
write(fstderr,*) "C6thl_Lscale0 = ", C6thl_Lscale0
write(fstderr,*) "C6rt_Lscale0 and C6thl_Lscale0 must be equal" &
// " when l_clip_semi_implicit is turned off" &
// " (default)."
err_code = clubb_fatal_error
endif ! C6rt_Lscale0 /= C6thl_Lscale0
endif ! .not. l_clip_semi_implicit
if ( C1 < zero ) then
write(fstderr,*) "C1 = ", C1
write(fstderr,*) "C1 must satisfy 0.0 <= C1"
err_code = clubb_fatal_error
end if
if ( C7 > one .or. C7 < zero ) then
write(fstderr,*) "C7 = ", C7
write(fstderr,*) "C7 must satisfy 0.0 <= C7 <= 1.0"
err_code = clubb_fatal_error
end if
if ( C7b > one .or. C7b < zero ) then
write(fstderr,*) "C7b = ", C7b
write(fstderr,*) "C7b must satisfy 0.0 <= C7b <= 1.0"
err_code = clubb_fatal_error
end if
if ( C11 > one .or. C11 < zero ) then
write(fstderr,*) "C11 = ", C11
write(fstderr,*) "C11 must satisfy 0.0 <= C11 <= 1.0"
err_code = clubb_fatal_error
end if
if ( C11b > one .or. C11b < zero ) then
write(fstderr,*) "C11b = ", C11b
write(fstderr,*) "C11b must satisfy 0.0 <= C11b <= 1.0"
err_code = clubb_fatal_error
end if
if ( C_wp2_splat < zero ) then
write(fstderr,*) "C_wp2_splat = ", C_wp2_splat
write(fstderr,*) "C_wp2_splat must satisfy C_wp2_splat >= 0"
err_code = clubb_fatal_error
end if
err_code_out = err_code
return
end subroutine setup_parameters
!=============================================================================
subroutine adj_low_res_nu &
( nzmax, grid_type, deltaz, & ! Intent(in)
momentum_heights, thermodynamic_heights, & ! Intent(in)
l_prescribed_avg_deltaz ) ! Intent(in)
! Description:
! Adjust the values of background eddy diffusivity based on
! vertical grid spacing.
! This code was made into a public subroutine so that it may be
! called multiple times per model run in scenarios where grid
! altitudes, and hence average grid spacing, change through space
! and/or time. This occurs, for example, when CLUBB is
! implemented in WRF. --ldgrant Jul 2010
!----------------------------------------------------------------------
use constants_clubb, only: &
fstderr ! Constant(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
implicit none
! Constant Parameters
! Flag for adjusting the values of the constant background eddy diffusivity
! coefficients based on the average vertical grid spacing. If this flag is
! turned off, the values of the various nu coefficients will remain as they
! are declared in the tunable_parameters.in file.
logical, parameter :: l_adj_low_res_nu = .true.
! The size of the average vertical grid spacing that serves as a threshold
! for when to increase the size of the background eddy diffusivity
! coefficients (nus) by a certain factor above what the background
! coefficients are specified to be in tunable_parameters.in. At any average
! grid spacing at or below this value, the values of the background
! diffusivities remain the same. However, at any average vertical grid
! spacing above this value, the values of the background eddy diffusivities
! are increased. Traditionally, the threshold grid spacing has been set to
! 40.0 meters. This is only relevant if l_adj_low_res_nu is turned on.
real( kind = core_rknd ), parameter :: &
grid_spacing_thresh = 40.0_core_rknd ! grid spacing threshold [m]
! Input Variables
! Grid definition
integer, intent(in) :: nzmax ! Vertical grid levels [#]
! If CLUBB is running on it's own, this option determines
! if it is using:
! 1) an evenly-spaced grid,
! 2) a stretched (unevenly-spaced) grid entered on the
! thermodynamic grid levels (with momentum levels set
! halfway between thermodynamic levels), or
! 3) a stretched (unevenly-spaced) grid entered on the
! momentum grid levels (with thermodynamic levels set
! halfway between momentum levels).
integer, intent(in) :: grid_type
real( kind = core_rknd ), intent(in) :: &
deltaz ! Change per height level [m]
! If the CLUBB parameterization is implemented in a host model,
! it needs to use the host model's momentum level altitudes
! and thermodynamic level altitudes.
! If the CLUBB model is running by itself, but is using a
! stretched grid entered on thermodynamic levels (grid_type = 2),
! it needs to use the thermodynamic level altitudes as input.
! If the CLUBB model is running by itself, but is using a
! stretched grid entered on momentum levels (grid_type = 3),
! it needs to use the momentum level altitudes as input.
real( kind = core_rknd ), intent(in), dimension(nzmax) :: &
momentum_heights, & ! Momentum level altitudes (input) [m]
thermodynamic_heights ! Thermodynamic level altitudes (input) [m]
logical, intent(in) :: &
l_prescribed_avg_deltaz ! used in adj_low_res_nu. If .true., avg_deltaz = deltaz
! Local Variables
real( kind = core_rknd ) :: avg_deltaz ! Average grid box height [m]
! The factor by which to multiply the coefficients of background eddy
! diffusivity if the grid spacing threshold is exceeded and l_adj_low_res_nu
! is turned on.
real( kind = core_rknd ),dimension(gr%nz) :: &
mult_factor_zt, & ! Uses gr%dzt for nu values on zt levels
mult_factor_zm ! Uses gr%dzm for nu values on zm levels
! Flag to enable nu values that are a function of grid spacing
logical, parameter :: l_nu_grid_dependent = .false.
integer :: k ! Loop variable
!--------------- Begin code -------------------------
if ( .not. allocated( nu1_vert_res_dep ) ) then
allocate( nu1_vert_res_dep(1:gr%nz) )
end if
if ( .not. allocated( nu2_vert_res_dep ) ) then
allocate( nu2_vert_res_dep(1:gr%nz) )
end if
if ( .not. allocated( nu6_vert_res_dep ) ) then
allocate( nu6_vert_res_dep(1:gr%nz) )
end if
if ( .not. allocated( nu8_vert_res_dep ) ) then
allocate( nu8_vert_res_dep(1:gr%nz) )
end if
if ( .not. allocated( nu9_vert_res_dep ) ) then
allocate( nu9_vert_res_dep(1:gr%nz) )
end if
if ( .not. allocated( nu10_vert_res_dep ) ) then
allocate( nu10_vert_res_dep(1:gr%nz) )
end if
if ( .not. allocated( nu_hm_vert_res_dep ) ) then
allocate( nu_hm_vert_res_dep(1:gr%nz) )
end if
! Flag for adjusting the values of the constant diffusivity coefficients
! based on the grid spacing. If this flag is turned off, the values of the
! various nu coefficients will remain as they are declared in the
! parameters.in file.
if ( l_adj_low_res_nu ) then
! ### Adjust Constant Diffusivity Coefficients Based On Grid Spacing ###
! All of the background coefficients of eddy diffusivity, as well as the
! constant coefficient for 4th-order hyper-diffusion, must be adjusted
! based on the size of the grid spacing. For a case that uses an
! evenly-spaced grid, the adjustment is based on the constant grid
! spacing deltaz. For a case that uses a stretched grid, the adjustment
! is based on avg_deltaz, which is the average grid spacing over the
! vertical domain.
if ( l_prescribed_avg_deltaz ) then
avg_deltaz = deltaz
else if ( grid_type == 3 ) then
! CLUBB is implemented in a host model, or is using grid_type = 3
! Find the average deltaz over the grid based on momentum level
! inputs.
avg_deltaz &
= ( momentum_heights(nzmax) - momentum_heights(1) ) &
/ real( nzmax - 1, kind = core_rknd )
else if ( grid_type == 1 ) then
! Evenly-spaced grid.
avg_deltaz = deltaz
else if ( grid_type == 2 ) then
! Stretched (unevenly-spaced) grid: stretched thermodynamic level
! input.
! Find the average deltaz over the stretched grid based on
! thermodynamic level inputs.
avg_deltaz &
= ( thermodynamic_heights(nzmax) - thermodynamic_heights(1) ) &
/ real( nzmax - 1, kind = core_rknd )
else
! Eric Raut added to remove compiler warning. (Obviously, this value is not used)
avg_deltaz = 0.0_core_rknd
write(fstderr,*) "Invalid grid_type:", grid_type
stop "Fatal error"
end if ! grid_type
! The nu's are chosen for deltaz <= 40 m. Looks like they must
! be adjusted for larger grid spacings (Vince Larson)
if( .not. l_nu_grid_dependent ) then
! Use a constant mult_factor so nu does not depend on grid spacing
if( avg_deltaz > grid_spacing_thresh ) then
mult_factor_zt = 1.0_core_rknd + mult_coef * log( avg_deltaz / grid_spacing_thresh )
mult_factor_zm = mult_factor_zt
else
mult_factor_zt = 1.0_core_rknd
mult_factor_zm = 1.0_core_rknd
end if
else ! l_nu_grid_dependent = .true.
! mult_factor will vary to create nu values that vary with grid spacing
do k = 1, gr%nz
if( gr%dzm(k) > grid_spacing_thresh ) then
mult_factor_zm(k) = 1.0_core_rknd + mult_coef * log( gr%dzm(k) / grid_spacing_thresh )
else
mult_factor_zm(k) = 1.0_core_rknd
end if
if( gr%dzt(k) > grid_spacing_thresh ) then
mult_factor_zt(k) = 1.0_core_rknd + mult_coef * log( gr%dzt(k) / grid_spacing_thresh )
else
mult_factor_zt(k) = 1.0_core_rknd
end if
end do
end if ! l_nu_grid_dependent
!mult_factor = 1.0_core_rknd + mult_coef * log( avg_deltaz / grid_spacing_thresh )
nu1_vert_res_dep = nu1 * mult_factor_zm
nu2_vert_res_dep = nu2 * mult_factor_zm
nu6_vert_res_dep = nu6 * mult_factor_zm
nu8_vert_res_dep = nu8 * mult_factor_zt
nu9_vert_res_dep = nu9 * mult_factor_zm
nu10_vert_res_dep = nu10 * mult_factor_zt !We're unsure of the grid
nu_hm_vert_res_dep = nu_hm * mult_factor_zt
else ! nu values are not adjusted
nu1_vert_res_dep = nu1
nu2_vert_res_dep = nu2
nu6_vert_res_dep = nu6
nu8_vert_res_dep = nu8
nu9_vert_res_dep = nu9
nu10_vert_res_dep = nu10
nu_hm_vert_res_dep = nu_hm
end if ! l_adj_low_res_nu
return
end subroutine adj_low_res_nu
#ifdef E3SM
!=============================================================================
subroutine clubb_param_readnl(filename)
! Description:
! Read clubb tunable parameters from namelist to override preset values
! To be called by clubb_readnl in clubb_intr.F90
! Author: Wuyin Lin
!-----------------------------------------------------------------------
use spmd_utils, only: masterproc
use namelist_utils, only: find_group_name
use units, only: getunit, freeunit
use cam_abortutils, only: endrun
use mpishorthand
implicit none
! Input variables
character(len=*), intent(in) :: filename
namelist /clubb_param_nl/ &
clubb_C1, &
clubb_C1b, &
clubb_C1c, &
clubb_C2rt, &
clubb_C2thl, &
clubb_C2rtthl, &
clubb_C4, &
clubb_C5, &
clubb_C6rt, &
clubb_C6rtb, &
clubb_C6rtc, &
clubb_C6thlb, &
clubb_C6thlc, &