-
Notifications
You must be signed in to change notification settings - Fork 0
/
advance_xm_wpxp_module.F90
6250 lines (5219 loc) · 274 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, &
wpxp_term_tp_lhs, &
wpxp_terms_ac_pr2_lhs, &
wpxp_term_pr1_lhs, &
wpxp_terms_bp_pr3_rhs, &
xm_correction_wpxp_cl, &
damp_coefficient, &
diagnose_upxp, &
error_prints_xm_wpxp
! 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
integer, parameter :: &
ndiags2 = 2, &
ndiags3 = 3, &
ndiags5 = 5
contains
!=============================================================================
subroutine advance_xm_wpxp( nz, ngrdcol, gr, dt, sigma_sqd_w, wm_zm, wm_zt, wp2, &
Lscale, wp3_on_wp2, wp3_on_wp2_zt, Kh_zt, Kh_zm, &
invrs_tau_C6_zm, tau_max_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, &
ice_supersat_frac, &
pdf_implicit_coefs_terms, &
um_forcing, vm_forcing, ug, vg, wpthvp, &
fcor, um_ref, vm_ref, up2, vp2, &
uprcp, vprcp, rc_coef, &
clubb_params, nu_vert_res_dep, &
iiPDF_type, &
penta_solve_method, &
tridiag_solve_method, &
l_predict_upwp_vpwp, &
l_diffuse_rtm_and_thlm, &
l_stability_correct_Kh_N2_zm, &
l_godunov_upwind_wpxp_ta, &
l_upwind_xm_ma, &
l_uv_nudge, &
l_tke_aniso, &
l_diag_Lscale_from_tau, &
l_use_C7_Richardson, &
l_brunt_vaisala_freq_moist, &
l_use_thvm_in_bv_freq, &
l_lmm_stepping, &
l_enable_relaxed_clipping, &
l_linearize_pbl_winds, &
l_mono_flux_lim_thlm, &
l_mono_flux_lim_rtm, &
l_mono_flux_lim_um, &
l_mono_flux_lim_vm, &
l_mono_flux_lim_spikefix, &
order_xm_wpxp, order_xp2_xpyp, order_wp2_wp3, &
stats_metadata, &
stats_zt, stats_zm, stats_sfc, &
rtm, wprtp, thlm, wpthlp, &
sclrm, wpsclrp, um, upwp, vm, vpwp, &
um_pert, vm_pert, upwp_pert, vpwp_pert )
! 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 parameter_indices, only: &
nparams, & ! Variable(s)
iC6rt, &
iC6rtb, &
iC6rtc, &
iC6thl, &
iC6thlb, &
iC6thlc, &
iC6rt_Lscale0, &
iC6thl_Lscale0, &
iC7, &
iC7b, &
iC7c, &
iC7_Lscale0, &
ic_K6, &
iwpxp_L_thresh, &
ialtitude_threshold, &
iC_uu_shr
use parameters_tunable, only: &
nu_vertical_res_dep ! Type(s)
use constants_clubb, only: &
fstderr, & ! Constant
one, &
one_half, &
zero, &
eps
use parameters_model, only: &
sclr_dim, & ! Variable(s)
ts_nudge
use grid_class, only: &
grid, & ! Type
ddzt ! Procedure(s)
use grid_class, only: &
zm2zt, & ! Procedure(s)
zt2zm
use model_flags, only: &
iiPDF_new, & ! Variable(s)
l_explicit_turbulent_adv_wpxp
use mono_flux_limiter, only: &
calc_turb_adv_range ! Procedure(s)
use pdf_parameter_module, only: &
implicit_coefs_terms ! Variable Type
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_metadata_type
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 stats_type, only: stats ! Type
implicit none
! -------------------- Input Variables --------------------
integer, intent(in) :: &
nz, &
ngrdcol
type (grid), target, intent(in) :: gr
real( kind = core_rknd ), intent(in) :: &
dt ! Timestep [s]
real( kind = core_rknd ), intent(in), dimension(ngrdcol,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
invrs_tau_C6_zm, & ! Inverse time-scale on mom. levels applied to C6 term [1/s]
tau_max_zm, & ! Max. allowable eddy dissipation time scale on m-levs [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(ngrdcol,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(ngrdcol,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 [-]
ice_supersat_frac
type(implicit_coefs_terms), 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(ngrdcol,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 ), dimension(ngrdcol,nz), intent(in) :: &
uprcp, & ! < u' r_c' > [(m kg)/(s kg)]
vprcp, & ! < v' r_c' > [(m kg)/(s kg)]
rc_coef ! Coefficient on X'r_c' in X'th_v' equation [K/(kg/kg)]
real( kind = core_rknd ), dimension(ngrdcol), intent(in) :: &
fcor ! Coriolis parameter [s^-1]
real( kind = core_rknd ), dimension(ngrdcol,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]
real( kind = core_rknd ), dimension(nparams), intent(in) :: &
clubb_params ! Array of CLUBB's tunable parameters [units vary]
type(nu_vertical_res_dep), intent(in) :: &
nu_vert_res_dep ! Vertical resolution dependent nu values
integer, intent(in) :: &
iiPDF_type, & ! Selected option for the two-component normal (double
! Gaussian) PDF type to use for the w, rt, and theta-l (or
! w, chi, and eta) portion of CLUBB's multivariate,
! two-component PDF.
penta_solve_method, & ! Method to solve then penta-diagonal system
tridiag_solve_method ! Specifier for method to solve tridiagonal systems
logical, intent(in) :: &
l_predict_upwp_vpwp, & ! Flag to predict <u'w'> and <v'w'> along with <u> and <v>
! alongside the advancement of <rt>, <w'rt'>, <thl>,
! <wpthlp>, <sclr>, and <w'sclr'> in subroutine
! advance_xm_wpxp. Otherwise, <u'w'> and <v'w'> are still
! approximated by eddy diffusivity when <u> and <v> are
! advanced in subroutine advance_windm_edsclrm.
l_diffuse_rtm_and_thlm, & ! This flag determines whether or not we want CLUBB to do
! diffusion on rtm and thlm
l_stability_correct_Kh_N2_zm, & ! This flag determines whether or not we want CLUBB to apply
! a stability correction
l_godunov_upwind_wpxp_ta, & ! This flag determines whether we want to use an upwind
! differencing approximation rather than a centered
! differencing for turbulent advection terms.
! It affects wpxp only.
l_upwind_xm_ma, & ! This flag determines whether we want to use an upwind
! differencing approximation rather than a centered
! differencing for turbulent or mean advection terms.
! It affects rtm, thlm, sclrm, um and vm.
l_uv_nudge, & ! For wind speed nudging
l_tke_aniso, & ! For anisotropic turbulent kinetic energy, i.e. TKE = 1/2
! (u'^2 + v'^2 + w'^2)
l_diag_Lscale_from_tau, & ! First diagnose dissipation time tau, and then diagnose the
! mixing length scale as Lscale = tau * tke
l_use_C7_Richardson, & ! Parameterize C7 based on Richardson number
l_brunt_vaisala_freq_moist, & ! Use a different formula for the Brunt-Vaisala frequency in
! saturated atmospheres (from Durran and Klemp, 1982)
l_use_thvm_in_bv_freq, & ! Use thvm in the calculation of Brunt-Vaisala frequency
l_lmm_stepping, & ! Apply Linear Multistep Method (LMM) Stepping
l_enable_relaxed_clipping, & ! Flag to relax clipping on wpxp in xm_wpxp_clipping_and_stats
l_linearize_pbl_winds, & ! Flag (used by E3SM) to linearize PBL winds
l_mono_flux_lim_thlm, & ! Flag to turn on monotonic flux limiter for thlm
l_mono_flux_lim_rtm, & ! Flag to turn on monotonic flux limiter for rtm
l_mono_flux_lim_um, & ! Flag to turn on monotonic flux limiter for um
l_mono_flux_lim_vm, & ! Flag to turn on monotonic flux limiter for vm
l_mono_flux_lim_spikefix ! Flag to implement monotonic flux limiter code that
! eliminates spurious drying tendencies at model top
integer, intent(in) :: &
order_xm_wpxp, &
order_xp2_xpyp, &
order_wp2_wp3
type (stats_metadata_type), intent(in) :: &
stats_metadata
! -------------------- Input/Output Variables --------------------
type (stats), target, dimension(ngrdcol), intent(inout) :: &
stats_zt, &
stats_zm, &
stats_sfc
real( kind = core_rknd ), intent(inout), dimension(ngrdcol,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]
real( kind = core_rknd ), intent(inout), dimension(ngrdcol,nz,sclr_dim) :: &
sclrm, & ! [Units vary]
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(ngrdcol,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]
! Variables used to track perturbed version of winds.
real( kind = core_rknd ), dimension(ngrdcol,nz), intent(inout) :: &
um_pert, & ! perturbed <u> [m/s]
vm_pert, & ! perturbed <v> [m/s]
upwp_pert, & ! perturbed <u'w'> [m^2/s^2]
vpwp_pert ! perturbed <v'w'> [m^2/s^2]
! -------------------- Local Variables --------------------
! Parameter Constants
logical, parameter :: &
l_iter = .true. ! True when the means and fluxes are prognosed
real( kind = core_rknd ) :: &
C6rt, & ! CLUBB tunable parameter C6rt
C6rtb, & ! CLUBB tunable parameter C6rtb
C6rtc, & ! CLUBB tunable parameter C6rtc
C6thl, & ! CLUBB tunable parameter C6thl
C6thlb, & ! CLUBB tunable parameter C6thlb
C6thlc, & ! CLUBB tunable parameter C6thlc
C6rt_Lscale0, & ! CLUBB tunable parameter C6rt_Lscale0
C6thl_Lscale0, & ! CLUBB tunable parameter C6thl_Lscale0
C7, & ! CLUBB tunable parameter C7
C7b, & ! CLUBB tunable parameter C7b
C7c, & ! CLUBB tunable parameter C7c
C7_Lscale0, & ! CLUBB tunable parameter C7_Lscale0
c_K6, & ! CLUBB tunable parameter c_K6
altitude_threshold, & ! CLUBB tunable parameter altitude_threshold
wpxp_L_thresh ! CLUBB tunable parameter wpxp_L_thresh
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
C6rt_Skw_fnc, C6thl_Skw_fnc, C7_Skw_fnc, C6_term
! Eddy Diffusion for wpthlp and wprtp.
real( kind = core_rknd ), dimension(ngrdcol,nz) :: Kw6 ! wpxp eddy diff. [m^2/s]
! 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(ngrdcol,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.
! Constant parameters as a function of Skw.
integer :: &
nrhs ! Number of RHS vectors
! Saved values of predictive fields, prior to being advanced, for use in
! print statements in case of fatal error.
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
rtm_old, & ! Saved value of r_t [kg/kg]
wprtp_old, & ! Saved value of w'r_t' [(kg/kg) m/s]
thlm_old, & ! Saved value of th_l [K]
wpthlp_old ! Saved value of w'th_l' [K m/s]
! Input/Output Variables
real( kind = core_rknd ), dimension(ngrdcol,nz,sclr_dim) :: &
sclrm_old, & ! Saved value of sclr [units vary]
wpsclrp_old ! Saved value of wpsclrp [units vary]
! Variables used to predict <u> and <u'w'>, as well as <v> and <v'w'>.
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
um_old, & ! Saved value of <u> [m/s]
upwp_old, & ! Saved value of <u'w'> [m^2/s^2]
vm_old, & ! Saved value of <v> [m/s]
vpwp_old ! Saved value of <v'w'> [m^2/s^2]
! LHS/RHS terms
real( kind = core_rknd ), dimension(ndiags3,ngrdcol,nz) :: &
lhs_diff_zm, & ! Diffusion term for w'x'
lhs_diff_zt, & ! Diffusion term for w'x'
lhs_ma_zt, & ! Mean advection contributions to lhs
lhs_ma_zm ! Mean advection contributions to lhs
real( kind = core_rknd ), dimension(ndiags3,ngrdcol,nz) :: &
lhs_ta_wprtp, & ! w'r_t' turbulent advection contributions to lhs
lhs_ta_wpthlp, & ! w'thl' turbulent advection contributions to lhs
lhs_ta_wpup, & ! w'u' turbulent advection contributions to lhs
lhs_ta_wpvp ! w'v' turbulent advection contributions to lhs
real( kind = core_rknd ), dimension(ndiags3,ngrdcol,nz,sclr_dim) :: &
lhs_ta_wpsclrp ! w'sclr' turbulent advection contributions to lhs
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
rhs_ta_wprtp, & ! w'r_t' turbulent advection contributions to rhs
rhs_ta_wpthlp, & ! w'thl' turbulent advection contributions to rhs
rhs_ta_wpup, & ! w'u' turbulent advection contributions to rhs
rhs_ta_wpvp ! w'v' turbulent advection contributions to rhs
real( kind = core_rknd ), dimension(ngrdcol,nz,sclr_dim) :: &
rhs_ta_wpsclrp ! w'sclr' turbulent advection contributions to rhs
real( kind = core_rknd ), dimension(ndiags2,ngrdcol,nz) :: &
lhs_tp, & ! Turbulent production terms of w'x'
lhs_ta_xm ! Turbulent advection terms of xm
real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
lhs_ac_pr2, & ! Accumulation of w'x' and w'x' pressure term 2
lhs_pr1_wprtp, & ! Pressure term 1 for w'r_t' for all grid levels
lhs_pr1_wpthlp, & ! Pressure term 1 for w'thl' for all grid levels
lhs_pr1_wpsclrp ! Pressure term 1 for w'sclr' for all grid levels
logical :: &
l_scalar_calc ! True if sclr_dim > 0
integer :: i, j, k
! Whether preturbed winds are being solved.
logical :: l_perturbed_wind
! -------------------- Begin Code --------------------
!$acc enter data create( C6rt_Skw_fnc, C6thl_Skw_fnc, C7_Skw_fnc, C6_term, Kw6, &
!$acc low_lev_effect, high_lev_effect, rtm_old, wprtp_old, thlm_old, &
!$acc wpthlp_old, um_old, upwp_old, vm_old, &
!$acc vpwp_old, lhs_diff_zm, lhs_diff_zt, lhs_ma_zt, lhs_ma_zm, &
!$acc lhs_ta_wprtp, lhs_ta_wpthlp, lhs_ta_wpup, lhs_ta_wpvp, &
!$acc rhs_ta_wprtp, rhs_ta_wpthlp, rhs_ta_wpup, &
!$acc rhs_ta_wpvp, lhs_tp, lhs_ta_xm, lhs_ac_pr2, &
!$acc lhs_pr1_wprtp, lhs_pr1_wpthlp )
!$acc enter data if( sclr_dim > 0 ) &
!$acc create( sclrm_old, wpsclrp_old, lhs_ta_wpsclrp, &
!$acc rhs_ta_wpsclrp, lhs_pr1_wpsclrp )
l_perturbed_wind = l_predict_upwp_vpwp .and. l_linearize_pbl_winds
! Check whether monotonic flux limiter flags are set appropriately
if ( clubb_at_least_debug_level( 0 ) ) then
if ( l_mono_flux_lim_rtm .and. .not. l_mono_flux_lim_spikefix ) then
write(fstderr,*) "l_mono_flux_lim_rtm=T with l_mono_flux_lim_spikefix=F can lead to spikes aloft."
err_code = clubb_fatal_error
return
end if
end if
! Check whether the passive scalars are present.
if ( sclr_dim > 0 ) then
l_scalar_calc = .true.
else
l_scalar_calc = .false.
end if
if ( iiPDF_type == iiPDF_new .and. ( .not. l_explicit_turbulent_adv_wpxp ) ) then
nrhs = 1
else
nrhs = 2 + sclr_dim
if ( l_predict_upwp_vpwp ) then
nrhs = nrhs + 2
if ( l_perturbed_wind ) then
nrhs = nrhs + 2
endif ! l_perturbed_wind
endif ! l_predict_upwp_vpwp
endif
! Save values of predictive fields to be printed in case of crash.
if ( l_lmm_stepping ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
rtm_old(i,k) = rtm(i,k)
wprtp_old(i,k) = wprtp(i,k)
thlm_old(i,k) = thlm(i,k)
wpthlp_old(i,k) = wpthlp(i,k)
end do
end do
!$acc end parallel loop
if ( sclr_dim > 0 ) then
!$acc parallel loop gang vector collapse(3) default(present)
do j = 1, sclr_dim
do k = 1, nz
do i = 1, ngrdcol
sclrm_old(i,k,j) = sclrm(i,k,j)
wpsclrp_old(i,k,j) = wpsclrp(i,k,j)
end do
end do
end do
!$acc end parallel loop
end if ! sclr_dim > 0
if ( l_predict_upwp_vpwp ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
um_old(i,k) = um(i,k)
upwp_old(i,k) = upwp(i,k)
vm_old(i,k) = vm(i,k)
vpwp_old(i,k) = vpwp(i,k)
end do
end do
!$acc end parallel loop
end if ! l_predict_upwp_vpwp
end if ! l_lmm_stepping
! Unpack CLUBB tunable parameters
C6rt = clubb_params(iC6rt)
C6thl = clubb_params(iC6thl)
altitude_threshold = clubb_params(ialtitude_threshold)
wpxp_L_thresh = clubb_params(iwpxp_L_thresh)
if ( .not. l_diag_Lscale_from_tau ) then
! Unpack CLUBB tunable parameters
C6rtb = clubb_params(iC6rtb)
C6rtc = clubb_params(iC6rtc)
C6thlb = clubb_params(iC6thlb)
C6thlc = clubb_params(iC6thlc)
C6rt_Lscale0 = clubb_params(iC6rt_Lscale0)
C6thl_Lscale0 = clubb_params(iC6thl_Lscale0)
! Compute C6 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
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C6rt_Skw_fnc(i,k) = C6rtb + ( C6rt - C6rtb ) &
* exp( -one_half * (Skw_zm(i,k)/C6rtc)**2 )
end do
end do
!$acc end parallel loop
else
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C6rt_Skw_fnc(i,k) = C6rtb
end do
end do
!$acc end parallel loop
end if
if ( abs(C6thl-C6thlb) > abs(C6thl+C6thlb)*eps/2 ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C6thl_Skw_fnc(i,k) = C6thlb + ( C6thl - C6thlb ) &
* exp( -one_half * (Skw_zm(i,k)/C6thlc)**2 )
end do
end do
!$acc end parallel loop
else
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C6thl_Skw_fnc(i,k) = C6thlb
end do
end do
!$acc end parallel loop
end if
! Damp C6 as a function of Lscale in stably stratified regions
call damp_coefficient( nz, ngrdcol, gr, C6rt, C6rt_Skw_fnc, &
C6rt_Lscale0, altitude_threshold, &
wpxp_L_thresh, Lscale, &
C6rt_Skw_fnc )
call damp_coefficient( nz, ngrdcol, gr, C6thl, C6thl_Skw_fnc, &
C6thl_Lscale0, altitude_threshold, &
wpxp_L_thresh, Lscale, &
C6thl_Skw_fnc )
else ! l_diag_Lscale_from_tau
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C6rt_Skw_fnc(i,k) = C6rt
C6thl_Skw_fnc(i,k) = C6thl
end do
end do
!$acc end parallel loop
endif ! .not. l_diag_Lscale_from_tau
! Compute C7_Skw_fnc
if ( l_use_C7_Richardson ) then
! New formulation based on Richardson number
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C7_Skw_fnc(i,k) = Cx_fnc_Richardson(i,k)
end do
end do
!$acc end parallel loop
else
! Unpack CLUBB tunable parameters
C7 = clubb_params(iC7)
C7b = clubb_params(iC7b)
C7c = clubb_params(iC7c)
C7_Lscale0 = clubb_params(iC7_Lscale0)
! Compute C7 as a function of Skw
if ( abs(C7-C7b) > abs(C7+C7b)*eps/2 ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C7_Skw_fnc(i,k) = C7b + ( C7 - C7b ) * exp( -one_half * (Skw_zm(i,k)/C7c)**2 )
end do
end do
!$acc end parallel loop
else
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C7_Skw_fnc(i,k) = C7b
end do
end do
!$acc end parallel loop
endif
! Damp C7 as a function of Lscale in stably stratified regions
call damp_coefficient( nz, ngrdcol, gr, C7, C7_Skw_fnc, &
C7_Lscale0, altitude_threshold, &
wpxp_L_thresh, Lscale, &
C7_Skw_fnc )
end if ! l_use_C7_Richardson
if ( stats_metadata%l_stats_samp ) then
!$acc update host( C7_Skw_fnc, C6rt_Skw_fnc, C6thl_Skw_fnc )
do i = 1, ngrdcol
call stat_update_var( stats_metadata%iC7_Skw_fnc, C7_Skw_fnc(i,:), & ! intent(in)
stats_zm(i) ) ! intent(inout)
call stat_update_var( stats_metadata%iC6rt_Skw_fnc, C6rt_Skw_fnc(i,:), & ! intent(in)
stats_zm(i) ) ! intent(inout
call stat_update_var( stats_metadata%iC6thl_Skw_fnc, C6thl_Skw_fnc(i,:), & ! intent(in)
stats_zm(i) ) ! intent(inout)
end do
end if
if ( clubb_at_least_debug_level( 0 ) ) then
! Assertion check for C7_Skw_fnc
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
if ( C7_Skw_fnc(i,k) > one .or. C7_Skw_fnc(i,k) < zero ) then
err_code = clubb_fatal_error
end if
end do
end do
!$acc end parallel loop
if ( err_code == clubb_fatal_error ) then
write(fstderr,*) "The C7_Skw_fnc variable is outside the valid range"
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
c_K6 = clubb_params(ic_K6)
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
Kw6(i,k) = c_K6 * Kh_zt(i,k)
end do
end do
!$acc end parallel loop
! 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( nz, ngrdcol, gr, dt, &
w_1_zm, w_2_zm, varnce_w_1_zm, varnce_w_2_zm, & ! intent(in)
mixt_frac_zm, & ! intent(in)
stats_metadata, & ! intent(in)
stats_zm, & ! intent(inout)
low_lev_effect, high_lev_effect ) ! intent(out)
! Calculate 1st pressure terms for w'r_t', w'thl', and w'sclr'.
call wpxp_term_pr1_lhs( nz, ngrdcol, C6rt_Skw_fnc, C6thl_Skw_fnc, C7_Skw_fnc, & ! Intent(in)
invrs_tau_C6_zm, l_scalar_calc, & ! Intent(in)
lhs_pr1_wprtp, lhs_pr1_wpthlp, lhs_pr1_wpsclrp ) ! Intent(out)
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
C6_term(i,k) = C6rt_Skw_fnc(i,k) * invrs_tau_C6_zm(i,k)
end do
end do
!$acc end parallel loop
if ( stats_metadata%l_stats_samp ) then
!$acc update host( C6_term )
do i = 1, ngrdcol
call stat_update_var( stats_metadata%iC6_term, C6_term(i,:), & ! intent(in)
stats_zm(i) ) ! intent(inout)
end do
end if
call calc_xm_wpxp_ta_terms( nz, ngrdcol, gr, wp2rtp, & ! intent(in)
wp2thlp, wp2sclrp, & ! intent(in)
rho_ds_zt, invrs_rho_ds_zm, rho_ds_zm, & ! intent(in)
sigma_sqd_w, wp3_on_wp2_zt, & ! intent(in)
pdf_implicit_coefs_terms, & ! intent(in)
iiPDF_type, & ! intent(in)
l_explicit_turbulent_adv_wpxp, l_predict_upwp_vpwp, & ! intent(in)
l_scalar_calc, & ! intent(in)
l_godunov_upwind_wpxp_ta, & ! intent(in)
stats_metadata, & ! intent(in)
stats_zt, & ! intent(inout)
lhs_ta_wprtp, lhs_ta_wpthlp, lhs_ta_wpup, & ! intent(out)
lhs_ta_wpvp, lhs_ta_wpsclrp, & ! intent(out)
rhs_ta_wprtp, rhs_ta_wpthlp, rhs_ta_wpup, & ! intent(out)
rhs_ta_wpvp, rhs_ta_wpsclrp ) ! intent(out)
! Calculate various terms that are the same between all LHS matricies
call calc_xm_wpxp_lhs_terms( nz, ngrdcol, gr, Kh_zm, wm_zm, wm_zt, wp2, & ! In
Kw6, C7_Skw_fnc, invrs_rho_ds_zt, & ! In
invrs_rho_ds_zm, rho_ds_zt, & ! In
rho_ds_zm, l_implemented, em, & ! In
Lscale, thlm, exner, rtm, rcm, p_in_Pa, thvm, & ! In
ice_supersat_frac, & ! In
clubb_params, nu_vert_res_dep, & ! In
l_diffuse_rtm_and_thlm, & ! In
l_stability_correct_Kh_N2_zm, & ! In
l_upwind_xm_ma, & ! In
l_brunt_vaisala_freq_moist, & ! In
l_use_thvm_in_bv_freq, & ! In
lhs_diff_zm, lhs_diff_zt, lhs_ma_zt, lhs_ma_zm, & ! Out
lhs_tp, lhs_ta_xm, lhs_ac_pr2 ) ! Out
! Setup and decompose matrix for each variable.
if ( ( iiPDF_type == iiPDF_new ) .and. ( .not. l_explicit_turbulent_adv_wpxp ) ) then
! LHS matrices are unique, multiple band solves required
call solve_xm_wpxp_with_multiple_lhs( nz, ngrdcol, gr, dt, l_iter, nrhs, wm_zt, wp2, & ! In
rtpthvp, rtm_forcing, wprtp_forcing, thlpthvp, & ! In
thlm_forcing, wpthlp_forcing, rho_ds_zm, & ! In
rho_ds_zt, invrs_rho_ds_zm, invrs_rho_ds_zt, & ! In
thv_ds_zm, rtp2, thlp2, l_implemented, & ! In
sclrpthvp, sclrm_forcing, sclrp2, & ! In
low_lev_effect, high_lev_effect, C7_Skw_fnc, & ! In
lhs_diff_zm, lhs_diff_zt, lhs_ma_zt, lhs_ma_zm, & ! In
lhs_ta_wprtp, lhs_ta_wpthlp, lhs_ta_wpsclrp, & ! In
rhs_ta_wprtp, rhs_ta_wpthlp, rhs_ta_wpsclrp, & ! In
lhs_tp, lhs_ta_xm, lhs_ac_pr2, lhs_pr1_wprtp, & ! In
lhs_pr1_wpthlp, lhs_pr1_wpsclrp, & ! In
penta_solve_method, & ! In
tridiag_solve_method, & ! In
l_predict_upwp_vpwp, & ! In
l_diffuse_rtm_and_thlm, & ! In
l_upwind_xm_ma, & ! In
l_tke_aniso, & ! In
l_enable_relaxed_clipping, & ! In
l_mono_flux_lim_thlm, & ! In
l_mono_flux_lim_rtm, & ! In
l_mono_flux_lim_um, & ! In
l_mono_flux_lim_vm, & ! In
l_mono_flux_lim_spikefix, & ! In
order_xm_wpxp, order_xp2_xpyp, order_wp2_wp3, & ! In
stats_metadata, & ! In
stats_zt, stats_zm, stats_sfc, & ! In
rtm, wprtp, thlm, wpthlp, sclrm, wpsclrp ) ! Out
else
! LHS matrices are equivalent, only one solve required
call solve_xm_wpxp_with_single_lhs( nz, ngrdcol, gr, dt, l_iter, nrhs, wm_zt, wp2, & ! In
invrs_tau_C6_zm, tau_max_zm, & ! In
rtpthvp, rtm_forcing, wprtp_forcing, thlpthvp, & ! In
thlm_forcing, wpthlp_forcing, rho_ds_zm, & ! In
rho_ds_zt, invrs_rho_ds_zm, invrs_rho_ds_zt, & ! In
thv_ds_zm, rtp2, thlp2, l_implemented, & ! In
sclrpthvp, sclrm_forcing, sclrp2, um_forcing, & ! In
vm_forcing, ug, vg, uprcp, vprcp, rc_coef, fcor, & ! In
up2, vp2, & ! In
low_lev_effect, high_lev_effect, & ! In
C6rt_Skw_fnc, C6thl_Skw_fnc, C7_Skw_fnc, & ! In
lhs_diff_zm, lhs_diff_zt, lhs_ma_zt, lhs_ma_zm, & ! In
lhs_ta_wprtp, & ! In
rhs_ta_wprtp, rhs_ta_wpthlp, rhs_ta_wpup, & ! In
rhs_ta_wpvp, rhs_ta_wpsclrp, & ! In
lhs_tp, lhs_ta_xm, lhs_ac_pr2, lhs_pr1_wprtp, & ! In
lhs_pr1_wpthlp, lhs_pr1_wpsclrp, & ! In
clubb_params(iC_uu_shr), & ! In
penta_solve_method, & ! In
tridiag_solve_method, & ! In
l_predict_upwp_vpwp, & ! In
l_diffuse_rtm_and_thlm, & ! In
l_upwind_xm_ma, & ! In
l_tke_aniso, & ! In
l_enable_relaxed_clipping, & ! In
l_perturbed_wind, & ! In
l_mono_flux_lim_thlm, & ! In
l_mono_flux_lim_rtm, & ! In
l_mono_flux_lim_um, & ! In
l_mono_flux_lim_vm, & ! In
l_mono_flux_lim_spikefix, & ! In
order_xm_wpxp, order_xp2_xpyp, order_wp2_wp3, & ! In
stats_metadata, & ! In
stats_zt, stats_zm, stats_sfc, & ! In
rtm, wprtp, thlm, wpthlp, & ! Out
sclrm, wpsclrp, um, upwp, vm, vpwp, & ! Out
um_pert, vm_pert, upwp_pert, vpwp_pert ) ! Out
end if ! ( ( iiPDF_type == iiPDF_new ) .and. ( .not. l_explicit_turbulent_adv_wpxp ) )
if ( l_lmm_stepping ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
thlm(i,k) = one_half * ( thlm_old(i,k) + thlm(i,k) )
rtm(i,k) = one_half * ( rtm_old(i,k) + rtm(i,k) )
wpthlp(i,k) = one_half * ( wpthlp_old(i,k) + wpthlp(i,k) )
wprtp(i,k) = one_half * ( wprtp_old(i,k) + wprtp(i,k) )
end do
end do
!$acc end parallel loop
if ( sclr_dim > 0 ) then
!$acc parallel loop gang vector collapse(3) default(present)
do j = 1, sclr_dim
do k = 1, nz
do i = 1, ngrdcol
sclrm(i,k,j) = one_half * ( sclrm_old(i,k,j) + sclrm(i,k,j) )
wpsclrp(i,k,j) = one_half * ( wpsclrp_old(i,k,j) + wpsclrp(i,k,j) )
end do
end do
end do
!$acc end parallel loop
endif ! sclr_dim > 0
if ( l_predict_upwp_vpwp ) then
!$acc parallel loop gang vector collapse(2) default(present)
do k = 1, nz
do i = 1, ngrdcol
um(i,k) = one_half * ( um_old(i,k) + um(i,k) )
vm(i,k) = one_half * ( vm_old(i,k) + vm(i,k) )
upwp(i,k) = one_half * ( upwp_old(i,k) + upwp(i,k) )
vpwp(i,k) = one_half * ( vpwp_old(i,k) + vpwp(i,k) )
end do
end do
!$acc end parallel loop
end if ! l_predict_upwp_vpwp
end if ! l_lmm_stepping
if ( clubb_at_least_debug_level( 0 ) ) then
if ( err_code == clubb_fatal_error ) then
!$acc update host( sigma_sqd_w, wm_zm, wm_zt, wp2, Lscale, wp3_on_wp2, &
!$acc wp3_on_wp2_zt, Kh_zt, Kh_zm, invrs_tau_C6_zm, Skw_zm, &
!$acc wp2rtp, rtpthvp, rtm_forcing, wprtp_forcing, rtm_ref, wp2thlp, &
!$acc thlpthvp, thlm_forcing, wpthlp_forcing, thlm_ref, rho_ds_zm, &
!$acc rho_ds_zt, invrs_rho_ds_zm, invrs_rho_ds_zt, thv_ds_zm, rtp2, &
!$acc thlp2, w_1_zm, w_2_zm, varnce_w_1_zm, varnce_w_2_zm, &
!$acc mixt_frac_zm, em, wp2sclrp, sclrpthvp, &
!$acc sclrm_forcing, sclrp2, exner, rcm, p_in_Pa, thvm, &
!$acc Cx_fnc_Richardson, um_forcing, vm_forcing, ug, vg, &
!$acc wpthvp, fcor, um_ref, vm_ref, up2, vp2, uprcp, vprcp, rc_coef, &
!$acc rtm, wprtp, thlm, wpthlp, sclrm, wpsclrp, um, upwp, vm, vpwp, &
!$acc rtm_old,wprtp_old, thlm_old, wpthlp_old, sclrm_old, wpsclrp_old, &
!$acc um_old, upwp_old, vm_old, vpwp_old )
do i = 1, ngrdcol
call error_prints_xm_wpxp( nz, gr%zm(i,:), gr%zt(i,:), & ! intent(in)
dt, sigma_sqd_w(i,:), wm_zm(i,:), wm_zt(i,:), wp2(i,:), & ! intent(in)
Lscale(i,:), wp3_on_wp2(i,:), wp3_on_wp2_zt(i,:), & ! intent(in)
Kh_zt(i,:), Kh_zm(i,:), invrs_tau_C6_zm(i,:), Skw_zm(i,:), & ! intent(in)
wp2rtp(i,:), rtpthvp(i,:), rtm_forcing(i,:), & ! intent(in)
wprtp_forcing(i,:), rtm_ref(i,:), wp2thlp(i,:), & ! intent(in)
thlpthvp(i,:), thlm_forcing(i,:), & ! intent(in)
wpthlp_forcing(i,:), thlm_ref(i,:), rho_ds_zm(i,:), & ! intent(in)
rho_ds_zt(i,:), invrs_rho_ds_zm(i,:), & ! intent(in)
invrs_rho_ds_zt(i,:), thv_ds_zm(i,:), rtp2(i,:), & ! intent(in)
thlp2(i,:), w_1_zm(i,:), w_2_zm(i,:), & ! intent(in)
varnce_w_1_zm(i,:), varnce_w_2_zm(i,:), & ! intent(in)
mixt_frac_zm(i,:), l_implemented, em(i,:), & ! intent(in)
wp2sclrp(i,:,:), sclrpthvp(i,:,:), sclrm_forcing(i,:,:), & ! intent(in)
sclrp2(i,:,:), exner(i,:), rcm(i,:), p_in_Pa(i,:), thvm(i,:), & ! intent(in)
Cx_fnc_Richardson(i,:), & ! intent(in)
pdf_implicit_coefs_terms, & ! intent(in)
um_forcing(i,:), vm_forcing(i,:), ug(i,:), vg(i,:), & ! intent(in)
wpthvp(i,:), fcor(i), um_ref(i,:), vm_ref(i,:), up2(i,:), & ! intent(in)
vp2(i,:), uprcp(i,:), vprcp(i,:), rc_coef(i,:), rtm(i,:), & ! intent(in)
wprtp(i,:), thlm(i,:), wpthlp(i,:), sclrm(i,:,:), wpsclrp(i,:,:), & ! intent(in)
um(i,:), upwp(i,:), vm(i,:), vpwp(i,:), rtm_old(i,:), & ! intent(in)
wprtp_old(i,:), thlm_old(i,:), wpthlp_old(i,:), & ! intent(in)
sclrm_old(i,:,:), wpsclrp_old(i,:,:), um_old(i,:), & ! intent(in)
upwp_old(i,:), vm_old(i,:), vpwp_old(i,:), & ! intent(in)
l_predict_upwp_vpwp, l_lmm_stepping ) ! intent(in)
end do
end if
end if
if ( rtm_sponge_damp_settings%l_sponge_damping ) then
!$acc update host( rtm, rtm_ref )
if ( stats_metadata%l_stats_samp ) then
do i = 1, ngrdcol
call stat_begin_update( nz, stats_metadata%irtm_sdmp, rtm(i,:) / dt, & ! intent(in)
stats_zt(i) ) ! intent(inout)
end do
end if
do i = 1, ngrdcol
rtm(i,:) = sponge_damp_xm( nz, dt, gr%zt(i,:), gr%zm(i,:), &
rtm_ref(i,:), rtm(i,:), rtm_sponge_damp_profile )
end do
if ( stats_metadata%l_stats_samp ) then
do i = 1, ngrdcol
call stat_end_update( nz, stats_metadata%irtm_sdmp, rtm(i,:) / dt, & ! intent(in)
stats_zt(i) ) ! intent(inout)
end do
end if
!$acc update device( rtm )
endif ! rtm_sponge_damp_settings%l_sponge_damping
if ( thlm_sponge_damp_settings%l_sponge_damping ) then
!$acc update host( thlm, thlm_ref )
if ( stats_metadata%l_stats_samp ) then