forked from ESCOMP/CLUBB_CESM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathadvance_xm_wpxp_module.F90
3975 lines (3202 loc) · 151 KB
/
advance_xm_wpxp_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
!-----------------------------------------------------------------------
! $Id$
!===============================================================================
module advance_xm_wpxp_module
! Description:
! Contains the CLUBB advance_xm_wpxp_module scheme.
! References:
! None
!-----------------------------------------------------------------------
implicit none
private ! Default scope
public :: advance_xm_wpxp
private :: xm_wpxp_lhs, &
xm_wpxp_rhs, &
xm_wpxp_solve, &
xm_wpxp_clipping_and_stats, &
xm_term_ta_lhs, &
xm_term_ta_lhs_all, &
wpxp_term_tp_lhs, &
wpxp_term_tp_lhs_all, &
wpxp_terms_ac_pr2_lhs, &
wpxp_terms_ac_pr2_lhs_all, &
wpxp_term_pr1_lhs, &
wpxp_term_pr1_lhs_all, &
wpxp_terms_bp_pr3_rhs, &
wpxp_terms_bp_pr3_rhs_all, &
xm_correction_wpxp_cl, &
damp_coefficient
! Parameter Constants
integer, parameter, private :: &
nsub = 2, & ! Number of subdiagonals in the LHS matrix
nsup = 2, & ! Number of superdiagonals in the LHS matrix
xm_wpxp_thlm = 1, & ! Named constant for thlm and wpthlp solving
xm_wpxp_rtm = 2, & ! Named constant for rtm and wprtp solving
xm_wpxp_scalar = 3, & ! Named constant for sclrm and wpsclrp solving
xm_wpxp_um = 4, & ! Named constant for optional um and upwp solving
xm_wpxp_vm = 5 ! Named constant for optional vm and vpwp solving
contains
!=============================================================================
subroutine advance_xm_wpxp( dt, sigma_sqd_w, wm_zm, wm_zt, wp2, &
Lscale, wp3_on_wp2, wp3_on_wp2_zt, Kh_zt, Kh_zm, &
tau_C6_zm, Skw_zm, wp2rtp, rtpthvp, rtm_forcing, &
wprtp_forcing, rtm_ref, wp2thlp, thlpthvp, &
thlm_forcing, wpthlp_forcing, thlm_ref, &
rho_ds_zm, rho_ds_zt, invrs_rho_ds_zm, &
invrs_rho_ds_zt, thv_ds_zm, rtp2, thlp2, &
w_1_zm, w_2_zm, varnce_w_1_zm, varnce_w_2_zm, &
mixt_frac_zm, l_implemented, em, wp2sclrp, &
sclrpthvp, sclrm_forcing, sclrp2, exner, rcm, &
p_in_Pa, thvm, Cx_fnc_Richardson, &
pdf_implicit_coefs_terms, &
um_forcing, vm_forcing, ug, vg, wpthvp, &
fcor, um_ref, vm_ref, up2, vp2, &
rtm, wprtp, thlm, wpthlp, &
sclrm, wpsclrp, um, upwp, vm, vpwp )
! Description:
! Advance the mean and flux terms by one timestep.
! References:
! https://arxiv.org/pdf/1711.03675v1.pdf#nameddest=url:wpxp_eqns
!
! Eqn. 16 & 17 on p. 3546 of
! ``A PDF-Based Model for Boundary Layer Clouds. Part I:
! Method and Model Description'' Golaz, et al. (2002)
! JAS, Vol. 59, pp. 3540--3551.
! See Also
! ``Equations for CLUBB'' Section 5:
! /Implicit solutions for the means and fluxes/
!-----------------------------------------------------------------------
use parameters_tunable, only: &
C6rt, & ! Variable(s)
C6rtb, &
C6rtc, &
C6thl, &
C6thlb, &
C6thlc, &
C7, &
C7b, &
C7c, &
c_K6, &
C6rt_Lscale0, &
C6thl_Lscale0, &
C7_Lscale0, &
wpxp_L_thresh
use constants_clubb, only: &
fstderr, & ! Constant
rt_tol, &
thl_tol, &
w_tol, &
w_tol_sqd, &
thl_tol_mfl, &
rt_tol_mfl, &
max_mag_correlation, &
one, &
one_half, &
zero, &
zero_threshold, &
eps
use parameters_model, only: &
sclr_dim, & ! Variable(s)
sclr_tol, &
ts_nudge
use grid_class, only: &
gr, & ! Variable(s)
ddzt ! Procedure(s)
use grid_class, only: &
zm2zt, & ! Procedure(s)
zt2zm
use model_flags, only: &
l_clip_semi_implicit, & ! Variable(s)
l_use_C7_Richardson, &
l_explicit_turbulent_adv_wpxp, &
l_upwind_wpxp_ta, &
l_predict_upwp_vpwp, &
l_uv_nudge
use mono_flux_limiter, only: &
calc_turb_adv_range ! Procedure(s)
use pdf_closure_module, only: &
iiPDF_new, & ! Variable(s)
iiPDF_ADG1, &
iiPDF_type
use pdf_parameter_module, only: &
implicit_coefs_terms ! Variable Type
use turbulent_adv_pdf, only: &
sgn_turbulent_velocity ! Procedure(s)
use clubb_precision, only: &
core_rknd ! Variable(s)
use error_code, only: &
clubb_at_least_debug_level, & ! Procedure
err_code, & ! Error Indicator
clubb_fatal_error ! Constants
use stats_type_utilities, only: &
stat_begin_update, & ! Procedure(s)
stat_end_update, &
stat_update_var
use stats_variables, only: &
stats_zt, &
stats_zm, &
irtm_matrix_condt_num, & ! Variables
ithlm_matrix_condt_num, &
irtm_sdmp, &
ithlm_sdmp, &
ium_sdmp, &
ivm_sdmp, &
ium_ndg, &
ivm_ndg, &
ium_ref, &
ivm_ref, &
ium_gf, &
ium_cf, &
ivm_gf, &
ivm_cf, &
ium_f, &
ivm_f, &
iupwp_pr4, &
ivpwp_pr4, &
iC7_Skw_fnc, &
iC6rt_Skw_fnc, &
iC6thl_Skw_fnc, &
icoef_wp2rtp_implicit, &
iterm_wp2rtp_explicit, &
icoef_wp2thlp_implicit, &
iterm_wp2thlp_explicit, &
l_stats_samp
use sponge_layer_damping, only: &
rtm_sponge_damp_settings, &
thlm_sponge_damp_settings, &
uv_sponge_damp_settings, &
rtm_sponge_damp_profile, &
thlm_sponge_damp_profile, &
uv_sponge_damp_profile, &
sponge_damp_xm ! Procedure(s)
use advance_helper_module, only: &
compute_Cx_fnc_Richardson ! Procedure
implicit none
! External
intrinsic :: exp, sqrt
! Parameter Constants
logical, parameter :: &
l_iter = .true. ! True when the means and fluxes are prognosed
! Input Variables
real( kind = core_rknd ), intent(in) :: &
dt ! Timestep [s]
real( kind = core_rknd ), intent(in), dimension(gr%nz) :: &
sigma_sqd_w, & ! sigma_sqd_w on momentum levels [-]
wm_zm, & ! w wind component on momentum levels [m/s]
wm_zt, & ! w wind component on thermodynamic levels [m/s]
wp2, & ! w'^2 (momentum levels) [m^2/s^2]
Lscale, & ! Turbulent mixing length [m]
em, & ! Turbulent Kinetic Energy (TKE) [m^2/s^2]
wp3_on_wp2, & ! Smoothed wp3 / wp2 on momentum levels [m/s]
wp3_on_wp2_zt, & ! Smoothed wp3 / wp2 on thermo. levels [m/s]
Kh_zt, & ! Eddy diffusivity on thermodynamic levels [m^2/s]
Kh_zm, & ! Eddy diffusivity on momentum levels
tau_C6_zm, & ! Time-scale tau on momentum levels applied to C6 term [s]
Skw_zm, & ! Skewness of w on momentum levels [-]
wp2rtp, & ! <w'^2 r_t'> (thermodynamic levels) [m^2/s^2 kg/kg]
rtpthvp, & ! r_t'th_v' (momentum levels) [(kg/kg) K]
rtm_forcing, & ! r_t forcing (thermodynamic levels) [(kg/kg)/s]
wprtp_forcing, & ! <w'r_t'> forcing (momentum levels) [(kg/kg)/s^2]
rtm_ref, & ! rtm for nudging [kg/kg]
wp2thlp, & ! <w'^2 th_l'> (thermodynamic levels) [m^2/s^2 K]
thlpthvp, & ! th_l'th_v' (momentum levels) [K^2]
thlm_forcing, & ! th_l forcing (thermodynamic levels) [K/s]
wpthlp_forcing, & ! <w'th_l'> forcing (momentum levels) [K/s^2]
thlm_ref, & ! thlm for nudging [K]
rho_ds_zm, & ! Dry, static density on momentum levels [kg/m^3]
rho_ds_zt, & ! Dry, static density on thermo. levels [kg/m^3]
invrs_rho_ds_zm, & ! Inv. dry, static density @ moment. levs. [m^3/kg]
invrs_rho_ds_zt, & ! Inv. dry, static density @ thermo. levs. [m^3/kg]
thv_ds_zm, & ! Dry, base-state theta_v on moment. levs. [K]
! Added for clipping by Vince Larson 29 Sep 2007
rtp2, & ! r_t'^2 (momentum levels) [(kg/kg)^2]
thlp2, & ! th_l'^2 (momentum levels) [K^2]
! End of Vince Larson's addition.
w_1_zm, & ! Mean w (1st PDF component) [m/s]
w_2_zm, & ! Mean w (2nd PDF component) [m/s]
varnce_w_1_zm, & ! Variance of w (1st PDF component) [m^2/s^2]
varnce_w_2_zm, & ! Variance of w (2nd PDF component) [m^2/s^2]
mixt_frac_zm ! Weight of 1st PDF component (Sk_w dependent) [-]
logical, intent(in) :: &
l_implemented ! Flag for CLUBB being implemented in a larger model.
! Additional variables for passive scalars
real( kind = core_rknd ), intent(in), dimension(gr%nz,sclr_dim) :: &
wp2sclrp, & ! <w'^2 sclr'> (thermodynamic levels) [Units vary]
sclrpthvp, & ! <sclr' th_v'> (momentum levels) [Units vary]
sclrm_forcing, & ! sclrm forcing (thermodynamic levels) [Units vary]
sclrp2 ! For clipping Vince Larson [Units vary]
real( kind = core_rknd ), intent(in), dimension(gr%nz) :: &
exner, & ! Exner function [-]
rcm, & ! cloud water mixing ratio, r_c [kg/kg]
p_in_Pa, & ! Air pressure [Pa]
thvm, & ! Virutal potential temperature [K]
Cx_fnc_Richardson ! Cx_fnc computed from Richardson_num [-]
type(implicit_coefs_terms), dimension(gr%nz), intent(in) :: &
pdf_implicit_coefs_terms ! Implicit coefs / explicit terms [units vary]
! Variables used to predict <u> and <u'w'>, as well as <v> and <v'w'>.
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
um_forcing, & ! <u> forcing term (thermodynamic levels) [m/s^2]
vm_forcing, & ! <v> forcing term (thermodynamic levels) [m/s^2]
ug, & ! <u> geostrophic wind (thermodynamic levels) [m/s]
vg, & ! <v> geostrophic wind (thermodynamic levels) [m/s]
wpthvp ! <w'thv'> (momentum levels) [m/s K]
real( kind = core_rknd ), intent(in) :: &
fcor ! Coriolis parameter [s^-1]
real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
um_ref, & ! Reference u wind component for nudging [m/s]
vm_ref, & ! Reference v wind component for nudging [m/s]
up2, & ! Variance of the u wind component [m^2/s^2]
vp2 ! Variance of the v wind component [m^2/s^2]
! Input/Output Variables
real( kind = core_rknd ), intent(inout), dimension(gr%nz) :: &
rtm, & ! r_t (total water mixing ratio) [kg/kg]
wprtp, & ! w'r_t' [(kg/kg) m/s]
thlm, & ! th_l (liquid water potential temperature) [K]
wpthlp ! w'th_l' [K m/s]
! Input/Output Variables
real( kind = core_rknd ), intent(inout), dimension(gr%nz,sclr_dim) :: &
sclrm, wpsclrp ! [Units vary]
! Variables used to predict <u> and <u'w'>, as well as <v> and <v'w'>.
real( kind = core_rknd ), intent(inout), dimension(gr%nz) :: &
um, & ! <u>: mean west-east horiz. velocity (thermo. levs.) [m/s]
upwp, & ! <u'w'>: momentum flux (momentum levels) [m^2/s^2]
vm, & ! <v>: mean south-north horiz. velocity (thermo. levs.) [m/s]
vpwp ! <v'w'>: momentum flux (momentum levels) [m^2/s^2]
! Local variables
real( kind = core_rknd ), dimension(nsup+nsub+1,2*gr%nz) :: &
lhs ! Implicit contributions to wpxp/xm (band diag. matrix) (LAPACK)
real( kind = core_rknd ), dimension(gr%nz) :: &
C6rt_Skw_fnc, C6thl_Skw_fnc, C7_Skw_fnc
! Additional variables for passive scalars
real( kind = core_rknd ), dimension(gr%nz,sclr_dim) :: &
wpsclrp_forcing ! <w'sclr'> forcing (momentum levels) [m/s{un vary}]
! Eddy Diffusion for wpthlp and wprtp.
real( kind = core_rknd ), dimension(gr%nz) :: Kw6 ! wpxp eddy diff. [m^2/s]
real( kind = core_rknd ), dimension(gr%nz) :: &
a1, & ! a_1 (momentum levels); See eqn. 24 in `Equations for CLUBB' [-]
a1_zt ! a_1 interpolated to thermodynamic levels [-]
! Variables for turbulent advection of predictive variances and covariances.
! <w'^2 rt'> = coef_wp2rtp_implicit * <w'rt'> + term_wp2rtp_explicit
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2rtp_implicit, & ! Coefficient that is multiplied by <w'rt'> [m/s]
term_wp2rtp_explicit ! Term that is on the RHS [m^2/s^2 kg/kg]
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2rtp_implicit_zm, & ! coef_wp2rtp_implicit interp. to m-levs. [m/s]
term_wp2rtp_explicit_zm ! term_wp2rtp_expl intrp m-levs [m^2/s^2 kg/kg]
! <w'^2 thl'> = coef_wp2thlp_implicit * <w'thl'> + term_wp2thlp_explicit
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2thlp_implicit, & ! Coef. that is multiplied by <w'thl'> [m/s]
term_wp2thlp_explicit ! Term that is on the RHS [m^2/s^2 K]
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2thlp_implicit_zm, & ! coef_wp2thlp_implicit interp. m-levs. [m/s]
term_wp2thlp_explicit_zm ! term_wp2thlp_expl interp. m-levs [m^2/s^2 K]
! <w'^2 sclr'> = coef_wp2sclrp_implicit * <w'sclr'> + term_wp2sclrp_explicit
real ( kind = core_rknd ), dimension(gr%nz,sclr_dim) :: &
coef_wp2sclrp_implicit, & ! Coef. that is multiplied by <w'sclr'> [m/s]
term_wp2sclrp_explicit ! Term that is on the RHS [m^2/s^2(un. vary)]
real ( kind = core_rknd ), dimension(gr%nz,sclr_dim) :: &
coef_wp2sclrp_implicit_zm, & ! coef_wp2sclrp_implicit interp. m-levs [m/s]
term_wp2sclrp_explicit_zm ! term_wp2sclrp_expl intrp zm [m^2/s^2(un v)]
! Sign of turbulent velocity (used for "upwind" turbulent advection)
real ( kind = core_rknd ), dimension(gr%nz) :: &
sgn_t_vel_wprtp, & ! Sign of the turbulent velocity for <w'rt'> [-]
sgn_t_vel_wpthlp ! Sign of the turbulent velocity for <w'thl'> [-]
real ( kind = core_rknd ), dimension(gr%nz,sclr_dim) :: &
sgn_t_vel_wpsclrp ! Sign of the turbulent velocity for <w'sclr'> [-]
! Variables used to predict <u> and <u'w'>, as well as <v> and <v'w'>.
! <w'^2 u'> = coef_wp2up_implicit * <u'w'> + term_wp2up_explicit
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2up_implicit, & ! Coefficient that is multiplied by <u'w'> [m/s]
term_wp2up_explicit ! Term that is on the RHS [m^3/s^3]
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2up_implicit_zm, & ! coef_wp2up_implicit interp. to m-levs. [m/s]
term_wp2up_explicit_zm ! term_wp2up_expl interp. to m-levs. [m^3/s^3]
! <w'^2 v'> = coef_wp2vp_implicit * <v'w'> + term_wp2vp_explicit
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2vp_implicit, & ! Coefficient that is multiplied by <v'w'> [m/s]
term_wp2vp_explicit ! Term that is on the RHS [m^3/s^3]
real ( kind = core_rknd ), dimension(gr%nz) :: &
coef_wp2vp_implicit_zm, & ! coef_wp2vp_implicit interp. to m-levs. [m/s]
term_wp2vp_explicit_zm ! term_wp2vp_expl interp. to m-levs. [m^3/s^3]
! Sign of turbulent velocity (used for "upwind" turbulent advection)
real ( kind = core_rknd ), dimension(gr%nz) :: &
sgn_t_vel_upwp, & ! Sign of the turbulent velocity for <u'w'> [-]
sgn_t_vel_vpwp ! Sign of the turbulent velocity for <v'w'> [-]
real( kind = core_rknd ), dimension(gr%nz) :: &
um_tndcy, & ! <u> forcing term + coriolis (t-levs) [m/s^2]
vm_tndcy, & ! <v> forcing term + coriolis (t-levs) [m/s^2]
upwp_forcing, & ! <u'w'> extra RHS pressure term (m-levs) [m^2/s^3]
vpwp_forcing, & ! <v'w'> extra RHS pressure term (m-levs) [m^2/s^3]
upthvp, & ! <u'thv'> (momentum levels) [m/s K]
vpthvp ! <v'thv'> (momentum levels) [m/s K]
! Variables used as part of the monotonic turbulent advection scheme.
! Find the lowermost and uppermost grid levels that can have an effect
! on the central thermodynamic level during the course of a time step,
! due to the effects of turbulent advection only.
integer, dimension(gr%nz) :: &
low_lev_effect, & ! Index of the lowest level that has an effect.
high_lev_effect ! Index of the highest level that has an effect.
! Variables used for clipping of w'x' due to correlation
! of w with x, such that:
! corr_(w,x) = w'x' / [ sqrt(w'^2) * sqrt(x'^2) ];
! -1 <= corr_(w,x) <= 1.
real( kind = core_rknd ), dimension(gr%nz) :: &
wpxp_upper_lim, & ! Keeps correlations from becoming greater than 1.
wpxp_lower_lim ! Keeps correlations from becoming less than -1.
real( kind = core_rknd ), dimension(gr%nz) :: dummy_1d ! Unreferenced array
real( kind = core_rknd ), allocatable, dimension(:,:) :: &
rhs, &! Right-hand sides of band diag. matrix. (LAPACK)
solution ! solution vectors of band diag. matrix. (LAPACK)
! Constant parameters as a function of Skw.
integer :: &
nrhs ! Number of RHS vectors
real( kind = core_rknd ) :: rcond
! Indices
integer :: i
!---------------------------------------------------------------------------
! ----- Begin Code -----
if ( l_clip_semi_implicit &
.or. ( ( iiPDF_type == iiPDF_new ) &
.and. ( .not. l_explicit_turbulent_adv_wpxp ) ) ) then
nrhs = 1
else
if ( l_predict_upwp_vpwp ) then
nrhs = 4+sclr_dim
else
nrhs = 2+sclr_dim
endif
endif
! Allocate rhs and solution vector
allocate( rhs(2*gr%nz,nrhs) )
allocate( solution(2*gr%nz,nrhs) )
! This is initialized solely for the purpose of avoiding a compiler
! warning about uninitialized variables.
dummy_1d = zero
! Compute C6 and C7 as a function of Skw
! The if...then is just here to save compute time
if ( abs(C6rt-C6rtb) > abs(C6rt+C6rtb)*eps/2 ) then
C6rt_Skw_fnc(1:gr%nz) = C6rtb + (C6rt-C6rtb) &
*EXP( -one_half * (Skw_zm(1:gr%nz)/C6rtc)**2 )
else
C6rt_Skw_fnc(1:gr%nz) = C6rtb
endif
if ( abs(C6thl-C6thlb) > abs(C6thl+C6thlb)*eps/2 ) then
C6thl_Skw_fnc(1:gr%nz) = C6thlb + (C6thl-C6thlb) &
*EXP( -one_half * (Skw_zm(1:gr%nz)/C6thlc)**2 )
else
C6thl_Skw_fnc(1:gr%nz) = C6thlb
endif
! Compute C7_Skw_fnc
if ( l_use_C7_Richardson ) then
! New formulation based on Richardson number
C7_Skw_fnc = Cx_fnc_Richardson
else
if ( abs(C7-C7b) > abs(C7+C7b)*eps/2 ) then
C7_Skw_fnc(1:gr%nz) = C7b + (C7-C7b) &
*EXP( -one_half * (Skw_zm(1:gr%nz)/C7c)**2 )
else
C7_Skw_fnc(1:gr%nz) = C7b
endif
! Damp C7 as a function of Lscale in stably stratified regions
C7_Skw_fnc = damp_coefficient( C7, C7_Skw_fnc, &
C7_Lscale0, wpxp_L_thresh, Lscale )
end if ! l_use_C7_Richardson
! Damp C6 as a function of Lscale in stably stratified regions
C6rt_Skw_fnc = damp_coefficient( C6rt, C6rt_Skw_fnc, &
C6rt_Lscale0, wpxp_L_thresh, Lscale )
C6thl_Skw_fnc = damp_coefficient( C6thl, C6thl_Skw_fnc, &
C6thl_Lscale0, wpxp_L_thresh, Lscale )
! C6rt_Skw_fnc = C6rt
! C6thl_Skw_fnc = C6thl
! C7_Skw_fnc = C7
if ( l_stats_samp ) then
call stat_update_var( iC7_Skw_fnc, C7_Skw_fnc, stats_zm )
call stat_update_var( iC6rt_Skw_fnc, C6rt_Skw_fnc, stats_zm )
call stat_update_var( iC6thl_Skw_fnc, C6thl_Skw_fnc, stats_zm )
end if
if ( clubb_at_least_debug_level( 0 ) ) then
! Assertion check for C7_Skw_fnc
if ( any( C7_Skw_fnc(:) > one ) .or. any( C7_Skw_fnc(:) < zero ) ) then
write(fstderr,*) "The C7_Skw_fnc variable is outside the valid range"
err_code = clubb_fatal_error
return
end if
end if
! Define the Coefficent of Eddy Diffusivity for the wpthlp and wprtp.
! Kw6 is used for wpthlp and wprtp, which are located on momentum levels.
! Kw6 is located on thermodynamic levels.
! Kw6 = c_K6 * Kh_zt
Kw6(1:gr%nz) = c_K6 * Kh_zt(1:gr%nz)
! Find the number of grid levels, both upwards and downwards, that can
! have an effect on the central thermodynamic level during the course of
! one time step due to turbulent advection. This is used as part of the
! monotonic turbulent advection scheme.
call calc_turb_adv_range( dt, w_1_zm, w_2_zm, varnce_w_1_zm, varnce_w_2_zm, & ! In
mixt_frac_zm, & ! In
low_lev_effect, high_lev_effect ) ! Out
! Set up the implicit coefficients and explicit terms for turbulent
! advection of <w'rt'>, <w'thl'>, and <w'sclr'>.
if ( l_explicit_turbulent_adv_wpxp ) then
! The turbulent advection of <w'x'> is handled explicitly.
! The <w'rt'> turbulent advection term is entirely explicit, as
! term_wp2rtp_explicit is equal to <w'^2 rt'> as calculated using PDF
! parameters, which is general for any PDF type. The value of
! <w'^2 rt'> is calculated on thermodynamic levels. The value of
! coef_wp2rtp_implicit is always 0.
coef_wp2rtp_implicit = zero
term_wp2rtp_explicit = wp2rtp
if ( l_upwind_wpxp_ta ) then
! Interpolate term_wp2rtp_explicit to momentum levels as
! term_wp2rtp_explicit_zm. The value of coef_wp2rtp_implicit_zm is
! always 0.
coef_wp2rtp_implicit_zm = zero
term_wp2rtp_explicit_zm = zt2zm( term_wp2rtp_explicit )
! Calculate the sign of the turbulent velocity for <w'rt'>.
sgn_t_vel_wprtp &
= sgn_turbulent_velocity( term_wp2rtp_explicit_zm, wprtp )
endif ! l_upwind_wpxp_ta
! The <w'thl'> turbulent advection term is entirely explicit, as
! term_wp2thlp_explicit is equal to <w'^2 thl'> as calculated using PDF
! parameters, which is general for any PDF type. The value of
! <w'^2 thl'> is calculated on thermodynamic levels. The value of
! coef_wp2thlp_implicit is always 0.
coef_wp2thlp_implicit = zero
term_wp2thlp_explicit = wp2thlp
if ( l_upwind_wpxp_ta ) then
! Interpolate term_wp2thlp_explicit to momentum levels as
! term_wp2thlp_explicit_zm. The value of coef_wp2thlp_implicit_zm is
! always 0.
coef_wp2thlp_implicit_zm = zero
term_wp2thlp_explicit_zm = zt2zm( term_wp2thlp_explicit )
! Calculate the sign of the turbulent velocity for <w'thl'>.
sgn_t_vel_wpthlp &
= sgn_turbulent_velocity( term_wp2thlp_explicit_zm, wpthlp )
endif ! l_upwind_wpxp_ta
do i = 1, sclr_dim, 1
! The <w'sclr'> turbulent advection term is entirely explicit, as
! term_wp2sclrp_explicit is equal to <w'^2 sclr'> as calculated
! using PDF parameters, which is general for any PDF type. The value
! of <w'^2 sclr'> is calculated on thermodynamic levels. The value of
! coef_wp2sclrp_implicit is always 0.
coef_wp2sclrp_implicit(:,i) = zero
term_wp2sclrp_explicit(:,i) = wp2sclrp(:,i)
if ( l_upwind_wpxp_ta ) then
! Interpolate term_wp2sclrp_explicit to momentum levels as
! term_wp2sclrp_explicit_zm. The value of
! coef_wp2sclrp_implicit_zm is always 0.
coef_wp2sclrp_implicit_zm(:,i) = zero
term_wp2sclrp_explicit_zm(:,i) &
= zt2zm( term_wp2sclrp_explicit(:,i) )
! Calculate the sign of the turbulent velocity for <w'sclr'>.
sgn_t_vel_wpsclrp(:,i) &
= sgn_turbulent_velocity( term_wp2sclrp_explicit_zm(:,i), &
wpsclrp(:,i) )
endif ! l_upwind_wpxp_ta
enddo ! i = 1, sclr_dim, 1
else ! .not. l_explicit_turbulent_adv_xpyp
! The turbulent advection of <w'x'> is handled implicitly or
! semi-implicitly.
if ( iiPDF_type == iiPDF_ADG1 ) then
! The ADG1 PDF is used.
! Calculate the implicit coefficients and explicit terms on
! thermodynamic grid levels.
! Calculate a_1.
! It is a variable that is a function of sigma_sqd_w (where
! sigma_sqd_w is located on momentum levels).
a1(1:gr%nz) = one / ( one - sigma_sqd_w(1:gr%nz) )
! Interpolate a_1 from momentum levels to thermodynamic levels. This
! will be used for the <w'x'> turbulent advection (ta) term.
a1_zt = max( zm2zt( a1 ), zero_threshold ) ! Positive def. quantity
! Implicit coefficient on <w'rt'> in <w'^2 rt'> equation.
coef_wp2rtp_implicit = a1_zt * wp3_on_wp2_zt
! The <w'rt'> turbulent advection term is entirely implicit, as
! <w'^2 rt'> = coef_wp2rtp_implicit * <w'rt'>. The value of
! term_wp2rtp_explicit is always 0.
term_wp2rtp_explicit = zero
if ( l_upwind_wpxp_ta ) then
! Calculate coef_wp2rtp_implicit on momentum levels as
! coef_wp2rtp_implicit_zm. The value of term_wp2rtp_explicit_zm is
! always 0.
coef_wp2rtp_implicit_zm = a1 * wp3_on_wp2
term_wp2rtp_explicit_zm = zero
! For ADG1, the sign of the turbulent velocity is the sign of
! <w'^3> / <w'^2>. For simplicity, the sign of turbulent velocity
! is set to wp3_on_wp2.
sgn_t_vel_wprtp = wp3_on_wp2
endif ! l_upwind_wpxp_ta
! Implicit coefficient on <w'thl'> in <w'^2 thl'> equation.
! For ADG1, this is the same as coef_wp2rtp_implicit.
coef_wp2thlp_implicit = coef_wp2rtp_implicit
! The <w'thl'> turbulent advection term is entirely implicit, as
! <w'^2 thl'> = coef_wp2thlp_implicit * <w'thl'>. The value of
! term_wp2thlp_explicit is always 0.
term_wp2thlp_explicit = zero
if ( l_upwind_wpxp_ta ) then
! Calculate coef_wp2thlp_implicit on momentum levels as
! coef_wp2thlp_implicit_zm. The value of term_wp2thlp_explicit_zm
! is always 0.
coef_wp2thlp_implicit_zm = coef_wp2rtp_implicit_zm
term_wp2thlp_explicit_zm = zero
! For ADG1, the sign of the turbulent velocity is the sign of
! <w'^3> / <w'^2>. For simplicity, the sign of turbulent velocity
! is set to wp3_on_wp2.
sgn_t_vel_wpthlp = wp3_on_wp2
endif ! l_upwind_wpxp_ta
do i = 1, sclr_dim, 1
! Implicit coefficient on <w'sclr'> in <w'^2 sclr'> equation.
! For ADG1, this is the same as coef_wp2rtp_implicit.
coef_wp2sclrp_implicit(:,i) = coef_wp2rtp_implicit
! The <w'sclr'> turbulent advection term is entirely implicit, as
! <w'^2 sclr'> = coef_wp2sclrp_implicit * <w'sclr'>. The value of
! term_wp2sclrp_explicit is always 0.
term_wp2sclrp_explicit(:,i) = zero
if ( l_upwind_wpxp_ta ) then
! Calculate coef_wp2sclrp_implicit on momentum levels as
! coef_wp2sclrp_implicit_zm. The value of
! term_wp2sclrp_explicit_zm is always 0.
coef_wp2sclrp_implicit_zm(:,i) = coef_wp2rtp_implicit_zm
term_wp2sclrp_explicit_zm(:,i) = zero
! For ADG1, the sign of the turbulent velocity is the sign of
! <w'^3> / <w'^2>. For simplicity, the sign of turbulent
! velocity is set to wp3_on_wp2.
sgn_t_vel_wpsclrp(:,i) = wp3_on_wp2
endif ! l_upwind_wpxp_ta
enddo ! i = 1, sclr_dim, 1
elseif ( iiPDF_type == iiPDF_new ) then
! The new PDF is used.
! Unpack the variables coef_wp2rtp_implicit, term_wp2rtp_explicit,
! coef_wp2thlp_implicit, and term_wp2thlp_explicit from
! pdf_implicit_coefs_terms. The PDF parameters and the resulting
! implicit coefficients and explicit terms are calculated on
! thermodynamic levels.
! Implicit coefficient on <w'rt'> in <w'^2 rt'> equation.
coef_wp2rtp_implicit = pdf_implicit_coefs_terms%coef_wp2rtp_implicit
! Explicit (RHS) term in <w'rt'> equation.
term_wp2rtp_explicit = pdf_implicit_coefs_terms%term_wp2rtp_explicit
if ( l_upwind_wpxp_ta ) then
! Interpolate coef_wp2rtp_implicit and term_wp2rtp_explicit
! to momentum levels as coef_wp2rtp_implicit_zm and
! term_wp2rtp_explicit_zm, respectively.
coef_wp2rtp_implicit_zm = zt2zm( coef_wp2rtp_implicit )
term_wp2rtp_explicit_zm = zt2zm( term_wp2rtp_explicit )
! Calculate the sign of the turbulent velocity for <w'rt'>.
sgn_t_vel_wprtp &
= sgn_turbulent_velocity( coef_wp2rtp_implicit_zm * wprtp &
+ term_wp2rtp_explicit_zm, wprtp )
endif ! l_upwind_wpxp_ta
! Implicit coefficient on <w'thl'> in <w'^2 thl'> equation.
coef_wp2thlp_implicit &
= pdf_implicit_coefs_terms%coef_wp2thlp_implicit
! Explicit (RHS) term in <w'thl'> equation.
term_wp2thlp_explicit &
= pdf_implicit_coefs_terms%term_wp2thlp_explicit
if ( l_upwind_wpxp_ta ) then
! Interpolate coef_wp2thlp_implicit and term_wp2thlp_explicit
! to momentum levels as coef_wp2thlp_implicit_zm and
! term_wp2thlp_explicit_zm, respectively.
coef_wp2thlp_implicit_zm = zt2zm( coef_wp2thlp_implicit )
term_wp2thlp_explicit_zm = zt2zm( term_wp2thlp_explicit )
! Calculate the sign of the turbulent velocity for <w'thl'>.
sgn_t_vel_wpthlp &
= sgn_turbulent_velocity( coef_wp2thlp_implicit_zm * wpthlp &
+ term_wp2thlp_explicit_zm, wpthlp )
endif ! l_upwind_wpxp_ta
do i = 1, sclr_dim, 1
! The code for the scalar variables will be set up later.
coef_wp2sclrp_implicit(:,i) = zero
term_wp2sclrp_explicit(:,i) = zero
if ( l_upwind_wpxp_ta ) then
! Interpolate coef_wp2sclrp_implicit and term_wp2sclrp_explicit
! to momentum levels as coef_wp2sclrp_implicit_zm and
! term_wp2sclrp_explicit_zm, respectively.
coef_wp2sclrp_implicit_zm(:,i) &
= zt2zm( coef_wp2sclrp_implicit(:,i) )
term_wp2sclrp_explicit_zm(:,i) &
= zt2zm( term_wp2sclrp_explicit(:,i) )
! Calculate the sign of the turbulent velocity for <w'sclr'>.
sgn_t_vel_wpsclrp(:,i) &
= sgn_turbulent_velocity( coef_wp2sclrp_implicit_zm(:,i) &
* wpsclrp(:,i) &
+ term_wp2sclrp_explicit_zm(:,i), &
wpsclrp(:,i) )
endif ! l_upwind_wpxp_ta
enddo ! i = 1, sclr_dim, 1
endif ! iiPDF_type
endif ! l_explicit_turbulent_adv_xpyp
if ( l_stats_samp ) then
call stat_update_var( icoef_wp2rtp_implicit, coef_wp2rtp_implicit, &
stats_zt )
call stat_update_var( iterm_wp2rtp_explicit, term_wp2rtp_explicit, &
stats_zt )
call stat_update_var( icoef_wp2thlp_implicit, coef_wp2thlp_implicit, &
stats_zt )
call stat_update_var( iterm_wp2thlp_explicit, term_wp2thlp_explicit, &
stats_zt )
endif
! Setup and decompose matrix for each variable.
if ( l_clip_semi_implicit &
.or. ( ( iiPDF_type == iiPDF_new ) &
.and. ( .not. l_explicit_turbulent_adv_wpxp ) ) ) then
! Compute the upper and lower limits of w'r_t' at every level,
! based on the correlation of w and r_t, such that:
! corr_(w,r_t) = w'r_t' / [ sqrt(w'^2) * sqrt(r_t'^2) ];
! -1 <= corr_(w,r_t) <= 1.
if ( l_clip_semi_implicit ) then
wpxp_upper_lim = max_mag_correlation * sqrt( wp2 * rtp2 )
wpxp_lower_lim = -wpxp_upper_lim
endif
! Compute the implicit portion of the r_t and w'r_t' equations.
! Build the left-hand side matrix.
call xm_wpxp_lhs( l_iter, dt, Kh_zm, wprtp, wm_zm, wm_zt, wp2, & ! In
coef_wp2rtp_implicit, coef_wp2rtp_implicit_zm, & ! In
sgn_t_vel_wprtp, Kw6, tau_C6_zm, C7_Skw_fnc, & ! In
C6rt_Skw_fnc, rho_ds_zm, rho_ds_zt, & ! In
invrs_rho_ds_zm, invrs_rho_ds_zt, & ! In
wpxp_upper_lim, wpxp_lower_lim, &! In
l_implemented, em, Lscale, thlm, exner, & ! In
rtm, rcm, p_in_Pa, thvm, & ! In
lhs ) ! Out
! Compute the explicit portion of the r_t and w'r_t' equations.
! Build the right-hand side vector.
call xm_wpxp_rhs( xm_wpxp_rtm, l_iter, dt, rtm, wprtp, & ! In
rtm_forcing, wprtp_forcing, C7_Skw_fnc, & ! In
rtpthvp, C6rt_Skw_fnc, tau_C6_zm, & ! In
coef_wp2rtp_implicit, coef_wp2rtp_implicit_zm, & ! In
term_wp2rtp_explicit, term_wp2rtp_explicit_zm, & ! In
sgn_t_vel_wprtp, rho_ds_zt, & ! In
rho_ds_zm, invrs_rho_ds_zm, thv_ds_zm, & ! In
wpxp_upper_lim, wpxp_lower_lim, & ! In
rhs(:,1) ) ! Out
! Solve r_t / w'r_t'
if ( l_stats_samp .and. irtm_matrix_condt_num > 0 ) then
call xm_wpxp_solve( nrhs, & ! Intent(in)
lhs, rhs, & ! Intent(inout)
solution, rcond ) ! Intent(out)
else
call xm_wpxp_solve( nrhs, & ! Intent(in)
lhs, rhs, & ! Intent(inout)
solution ) ! Intent(out)
endif
if ( clubb_at_least_debug_level( 0 ) ) then
if ( err_code == clubb_fatal_error ) then
write(fstderr,'(a)') "Mean total water & total water flux LU decomp. failed"
return
end if
end if
call xm_wpxp_clipping_and_stats &
( xm_wpxp_rtm, dt, wp2, rtp2, wm_zt, & ! Intent(in)
rtm_forcing, rho_ds_zm, rho_ds_zt, & ! Intent(in)
invrs_rho_ds_zm, invrs_rho_ds_zt, & ! Intent(in)
rt_tol**2, rt_tol, rcond, & ! Intent(in)
low_lev_effect, high_lev_effect, & ! Intent(in)
l_implemented, solution(:,1), & ! Intent(in)
rtm, rt_tol_mfl, wprtp ) ! Intent(inout)
if ( clubb_at_least_debug_level( 0 ) ) then
if ( err_code == clubb_fatal_error ) then
write(fstderr,'(a)') "rtm monotonic flux limiter: tridag failed"
return
end if
end if
! Compute the upper and lower limits of w'th_l' at every level,
! based on the correlation of w and th_l, such that:
! corr_(w,th_l) = w'th_l' / [ sqrt(w'^2) * sqrt(th_l'^2) ];
! -1 <= corr_(w,th_l) <= 1.
if ( l_clip_semi_implicit ) then
wpxp_upper_lim = max_mag_correlation * sqrt( wp2 * thlp2 )
wpxp_lower_lim = -wpxp_upper_lim
endif
! Compute the implicit portion of the th_l and w'th_l' equations.
! Build the left-hand side matrix.
call xm_wpxp_lhs( l_iter, dt, Kh_zm, wpthlp, wm_zm, wm_zt, wp2, & ! In
coef_wp2thlp_implicit, coef_wp2thlp_implicit_zm, & ! In
sgn_t_vel_wpthlp, Kw6, tau_C6_zm, C7_Skw_fnc, & ! In
C6thl_Skw_fnc, rho_ds_zm, rho_ds_zt, & ! In
invrs_rho_ds_zm, invrs_rho_ds_zt, & ! In
wpxp_upper_lim, wpxp_lower_lim, &! In
l_implemented, em, Lscale, thlm, exner, & ! In
rtm, rcm, p_in_Pa, thvm, & ! In
lhs ) ! Out
! Compute the explicit portion of the th_l and w'th_l' equations.
! Build the right-hand side vector.
call xm_wpxp_rhs( xm_wpxp_thlm, l_iter, dt, thlm, wpthlp, & ! In
thlm_forcing, wpthlp_forcing, C7_Skw_fnc, & ! In
thlpthvp, C6thl_Skw_fnc, tau_C6_zm, & ! In
coef_wp2thlp_implicit, coef_wp2thlp_implicit_zm, & ! In
term_wp2thlp_explicit, term_wp2thlp_explicit_zm, & ! In
sgn_t_vel_wpthlp, rho_ds_zt, & ! In
rho_ds_zm, invrs_rho_ds_zm, thv_ds_zm, & ! In
wpxp_upper_lim, wpxp_lower_lim, & ! In
rhs(:,1) ) ! Out
! Solve for th_l / w'th_l'
if ( l_stats_samp .and. ithlm_matrix_condt_num > 0 ) then
call xm_wpxp_solve( nrhs, & ! Intent(in)
lhs, rhs, & ! Intent(inout)
solution, rcond ) ! Intent(out)
else
call xm_wpxp_solve( nrhs, & ! Intent(in)
lhs, rhs, & ! Intent(inout)
solution ) ! Intent(out)
endif
if ( clubb_at_least_debug_level( 0 ) ) then
if ( err_code == clubb_fatal_error ) then
write(fstderr,'(a)') "Liquid pot. temp & thetal flux LU decomp. failed"
return
end if
end if
call xm_wpxp_clipping_and_stats &
( xm_wpxp_thlm, dt, wp2, thlp2, wm_zt, & ! Intent(in)
thlm_forcing, rho_ds_zm, rho_ds_zt, & ! Intent(in)
invrs_rho_ds_zm, invrs_rho_ds_zt, & ! Intent(in)
thl_tol**2, thl_tol, rcond, & ! Intent(in)
low_lev_effect, high_lev_effect, & ! Intent(in)
l_implemented, solution(:,1), & ! Intent(in)
thlm, thl_tol_mfl, wpthlp ) ! Intent(inout)
if ( clubb_at_least_debug_level( 0 ) ) then
if ( err_code == clubb_fatal_error ) then
write(fstderr,'(a)') "thlm monotonic flux limiter: tridag failed"
return
end if
end if
! Solve sclrm / wpsclrp
! If sclr_dim is 0, then this loop will execute 0 times.
! ---> h1g, 2010-06-15
! scalar transport, e.g, droplet and ice number concentration
! are handled in " advance_sclrm_Nd_module.F90 "
#ifdef GFDL
do i = 1, 0, 1
#else
do i = 1, sclr_dim, 1
#endif
! <--- h1g, 2010-06-15
! Compute the upper and lower limits of w'sclr' at every level,
! based on the correlation of w and sclr, such that:
! corr_(w,sclr) = w'sclr' / [ sqrt(w'^2) * sqrt(sclr'^2) ];
! -1 <= corr_(w,sclr) <= 1.
if ( l_clip_semi_implicit ) then
wpxp_upper_lim(:) = max_mag_correlation * sqrt( wp2(:) * sclrp2(:,i) )
wpxp_lower_lim(:) = -wpxp_upper_lim(:)
endif
! Set <w'sclr'> forcing to 0 unless unless testing the wpsclrp code
! using wprtp or wpthlp (then use wprtp_forcing or wpthlp_forcing).
wpsclrp_forcing(:,i) = zero
! Compute the implicit portion of the sclr and w'sclr' equations.
! Build the left-hand side matrix.
call xm_wpxp_lhs( l_iter, dt, Kh_zm, wpsclrp(:,i), & ! In
wm_zm, wm_zt, wp2, & ! In
coef_wp2sclrp_implicit(:,i), & ! In
coef_wp2sclrp_implicit_zm(:,i), & ! In
sgn_t_vel_wpsclrp(:,i), Kw6, tau_C6_zm, & ! In
C7_Skw_fnc, C6rt_Skw_fnc, rho_ds_zm, rho_ds_zt, & ! In
invrs_rho_ds_zm, invrs_rho_ds_zt, & ! In
wpxp_upper_lim, wpxp_lower_lim, &! In
l_implemented, em, Lscale, thlm, exner, & ! In
rtm, rcm, p_in_Pa, thvm, & ! In
lhs ) ! Out
! Compute the explicit portion of the sclrm and w'sclr' equations.
! Build the right-hand side vector.
call xm_wpxp_rhs( xm_wpxp_scalar, l_iter, dt, sclrm(:,i), & ! In
wpsclrp(:,i), sclrm_forcing(:,i), &
wpsclrp_forcing(:,i), C7_Skw_fnc, & ! In
sclrpthvp(:,i), C6rt_Skw_fnc, tau_C6_zm, & ! In
coef_wp2sclrp_implicit(:,i), & ! In
coef_wp2sclrp_implicit_zm(:,i), & ! In
term_wp2sclrp_explicit(:,i), & ! In
term_wp2sclrp_explicit_zm(:,i), & ! In
sgn_t_vel_wpsclrp(:,i), rho_ds_zt, & ! In
rho_ds_zm, invrs_rho_ds_zm, thv_ds_zm, & ! In
wpxp_upper_lim, wpxp_lower_lim, & ! In
rhs(:,1) ) ! Out
! Solve for sclrm / w'sclr'
call xm_wpxp_solve( nrhs, & ! Intent(in)
lhs, rhs, & ! Intent(inout)
solution ) ! Intent(out)
if ( clubb_at_least_debug_level( 0 ) ) then
if ( err_code == clubb_fatal_error ) then
write(fstderr,*) "Passive scalar # ", i, " LU decomp. failed."
return