Split cloud components according to saturation pressures; source theta from latent heat.
24 int SAM_moisture_type = 1;
27 SAM_moisture_type = 2;
30 auto domain =
m_geom.Domain();
31 int i_lo = domain.smallEnd(0);
32 int i_hi = domain.bigEnd(0);
33 int j_lo = domain.smallEnd(1);
34 int j_hi = domain.bigEnd(1);
48 auto tbx = mfi.tilebox();
49 if (tbx.smallEnd(0) == i_lo) { tbx.growLo(0,-
m_real_width); }
50 if (tbx.bigEnd(0) == i_hi) { tbx.growHi(0,-
m_real_width); }
51 if (tbx.smallEnd(1) == j_lo) { tbx.growLo(1,-
m_real_width); }
52 if (tbx.bigEnd(1) == j_hi) { tbx.growHi(1,-
m_real_width); }
54 ParallelFor(tbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k)
63 Real delta_qv, delta_qc, delta_qi;
71 if (SAM_moisture_type == 1){
73 if (tabs_array(i,j,k) >=
tbgmax) {
75 delta_qi = qci_array(i,j,k);
76 qci_array(i,j,k) = 0.0;
77 qcl_array(i,j,k) += delta_qi;
78 tabs_array(i,j,k) -= fac_fus * delta_qi;
79 pres_array(i,j,k) = rho_array(i,j,k) *
R_d * tabs_array(i,j,k)
80 * (1.0 +
R_v/
R_d * qv_array(i,j,k));
81 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), pres_array(i,j,k), rdOcp);
82 pres_array(i,j,k) *= 0.01;
85 else if (tabs_array(i,j,k) <=
tbgmin) {
87 delta_qc = qcl_array(i,j,k);
88 qcl_array(i,j,k) = 0.0;
89 qci_array(i,j,k) += delta_qc;
90 tabs_array(i,j,k) += fac_fus * delta_qc;
91 pres_array(i,j,k) = rho_array(i,j,k) *
R_d * tabs_array(i,j,k)
92 * (1.0 +
R_v/
R_d * qv_array(i,j,k));
93 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), pres_array(i,j,k), rdOcp);
94 pres_array(i,j,k) *= 0.01;
98 omn = an*tabs_array(i,j,k)-bn;
99 delta_qc = qcl_array(i,j,k) - qn_array(i,j,k) * omn;
100 delta_qi = qci_array(i,j,k) - qn_array(i,j,k) * (1.0 - omn);
101 qcl_array(i,j,k) = qn_array(i,j,k) * omn;
102 qci_array(i,j,k) = qn_array(i,j,k) * (1.0 - omn);
103 tabs_array(i,j,k) += fac_fus * delta_qc;
104 pres_array(i,j,k) = rho_array(i,j,k) *
R_d * tabs_array(i,j,k)
105 * (1.0 +
R_v/
R_d * qv_array(i,j,k));
106 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), pres_array(i,j,k), rdOcp);
107 pres_array(i,j,k) *= 0.01;
110 else if (SAM_moisture_type == 2)
113 delta_qc = qcl_array(i,j,k) - qn_array(i,j,k);
115 qcl_array(i,j,k) = qn_array(i,j,k);
116 qci_array(i,j,k) = 0.0;
117 tabs_array(i,j,k) += fac_cond * delta_qc;
118 pres_array(i,j,k) = rho_array(i,j,k) *
R_d * tabs_array(i,j,k)
119 * (1.0 +
R_v/
R_d * qv_array(i,j,k));
120 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), pres_array(i,j,k), rdOcp);
121 pres_array(i,j,k) *= 0.01;
125 erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsatw);
126 erf_qsati(tabs_array(i,j,k), pres_array(i,j,k), qsati);
127 qsat = omn * qsatw + (1.0-omn) * qsati;
130 if (qt_array(i,j,k) > qsat) {
133 tabs_array(i,j,k) =
NewtonIterSat(i, j, k , SAM_moisture_type ,
134 fac_cond , fac_fus , fac_sub ,
136 tabs_array, pres_array,
137 qv_array , qcl_array , qci_array,
138 qn_array , qt_array);
141 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), 100.0*pres_array(i,j,k), rdOcp);
152 delta_qv = qcl_array(i,j,k) + qci_array(i,j,k);
153 delta_qc = qcl_array(i,j,k);
154 delta_qi = qci_array(i,j,k);
157 qv_array(i,j,k) += delta_qv;
158 qcl_array(i,j,k) = 0.0;
159 qci_array(i,j,k) = 0.0;
160 qn_array(i,j,k) = 0.0;
161 qt_array(i,j,k) = qv_array(i,j,k);
164 tabs_array(i,j,k) -= fac_cond * delta_qc + fac_sub * delta_qi;
167 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), 100.0*pres_array(i,j,k), rdOcp);
170 erf_qsatw(tabs_array(i,j,k), pres_array(i,j,k), qsatw);
171 erf_qsati(tabs_array(i,j,k), pres_array(i,j,k), qsati);
172 qsat = omn * qsatw + (1.0-omn) * qsati;
173 if (qt_array(i,j,k) > qsat) {
176 tabs_array(i,j,k) =
NewtonIterSat(i, j, k , SAM_moisture_type ,
177 fac_cond , fac_fus , fac_sub ,
179 tabs_array, pres_array,
180 qv_array , qcl_array , qci_array,
181 qn_array , qt_array);
184 theta_array(i,j,k) =
getThgivenTandP(tabs_array(i,j,k), 100.0*pres_array(i,j,k), rdOcp);
constexpr amrex::Real R_v
Definition: ERF_Constants.H:11
constexpr amrex::Real tbgmax
Definition: ERF_Constants.H:32
constexpr amrex::Real tbgmin
Definition: ERF_Constants.H:31
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getThgivenTandP(const amrex::Real T, const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:18
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:166
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsati(amrex::Real t, amrex::Real p, amrex::Real &qsati)
Definition: ERF_MicrophysicsUtils.H:155
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Array< FabPtr, MicVar::NumVars > mic_fab_vars
Definition: ERF_SAM.H:329
amrex::Real m_rdOcp
Definition: ERF_SAM.H:318
amrex::Real m_fac_fus
Definition: ERF_SAM.H:316
bool m_do_cond
Definition: ERF_SAM.H:319
amrex::Real m_fac_cond
Definition: ERF_SAM.H:315
amrex::Geometry m_geom
Definition: ERF_SAM.H:300
amrex::Real m_fac_sub
Definition: ERF_SAM.H:317
int m_real_width
Definition: ERF_SAM.H:303
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real NewtonIterSat(int &i, int &j, int &k, const int &SAM_moisture_type, const amrex::Real &fac_cond, const amrex::Real &, const amrex::Real &fac_sub, const amrex::Real &an, const amrex::Real &bn, const amrex::Array4< amrex::Real > &tabs_array, const amrex::Array4< amrex::Real > &pres_array, const amrex::Array4< amrex::Real > &qv_array, const amrex::Array4< amrex::Real > &qc_array, const amrex::Array4< amrex::Real > &qi_array, const amrex::Array4< amrex::Real > &qn_array, const amrex::Array4< amrex::Real > &qt_array)
Definition: ERF_SAM.H:160
@ pres
Definition: ERF_SAM.H:33
@ qci
Definition: ERF_SAM.H:39
@ qv
Definition: ERF_SAM.H:37
@ rho
Definition: ERF_SAM.H:30
@ qt
Definition: ERF_SAM.H:35
@ qn
Definition: ERF_SAM.H:36
@ theta
Definition: ERF_SAM.H:31
@ qcl
Definition: ERF_SAM.H:38
@ tabs
Definition: ERF_SAM.H:32
MoistureType moisture_type
Definition: ERF_DataStruct.H:1171