54 #ifndef ERF_RICHARDSON_NUMBER_H_
55 #define ERF_RICHARDSON_NUMBER_H_
57 #include <AMReX_GpuQualifiers.H>
58 #include <AMReX_REAL.H>
81 const amrex::Array4<const amrex::Real>& cell_data,
129 const amrex::Array4<const amrex::Real>& cell_data,
134 const int qv_index = moisture_indices.
qv;
135 const int qc_index = moisture_indices.
qc;
139 bool saturated =
IsSaturated(i, j, k, cell_data, qc_index);
144 if (saturated && qc_index >= 0) {
149 inv_theta = 1.0 / theta_l;
150 dtheta_dz = 0.5 * (theta_l_kp1 - theta_l_km1) * dzInv;
152 }
else if (qv_index >= 0) {
157 inv_theta = 1.0 / theta_v;
158 dtheta_dz = 0.5 * (theta_v_kp1 - theta_v_km1) * dzInv;
162 amrex::Real theta_kp1 = cell_data(i, j, k+1, rhoTheta_comp) / cell_data(i, j, k+1, rho_comp);
163 amrex::Real theta_km1 = cell_data(i, j, k-1, rhoTheta_comp) / cell_data(i, j, k-1, rho_comp);
165 inv_theta = 1.0 /
theta;
166 dtheta_dz = 0.5 * (theta_kp1 - theta_km1) * dzInv;
169 return const_grav * inv_theta * dtheta_dz;
203 const amrex::Array4<const amrex::Real>& u,
204 const amrex::Array4<const amrex::Real>& v)
206 amrex::Real u_c_km1 = 0.5 * (u(i, j, k-1) + u(i+1, j, k-1));
207 amrex::Real u_c_kp1 = 0.5 * (u(i, j, k+1) + u(i+1, j, k+1));
209 amrex::Real v_c_km1 = 0.5 * (v(i, j, k-1) + v(i, j+1, k-1));
210 amrex::Real v_c_kp1 = 0.5 * (v(i, j, k+1) + v(i, j+1, k+1));
212 amrex::Real du_dz = 0.5 * (u_c_kp1 - u_c_km1) * dzInv;
213 amrex::Real dv_dz = 0.5 * (v_c_kp1 - v_c_km1) * dzInv;
215 amrex::Real S2_vert = du_dz * du_dz + dv_dz * dv_dz;
232 amrex::Real S2_safe = amrex::max(S2_vert, 1.0e-10);
233 return N2_moist / S2_safe;
250 }
else if (Ri > 0.0) {
251 return std::sqrt(1.0 - Ri / Ri_crit);
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real GetThetal(int i, int j, int k, const amrex::Array4< amrex::Real const > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_MoistUtils.H:101
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real GetThetav(int i, int j, int k, const amrex::Array4< amrex::Real const > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_MoistUtils.H:72
AMREX_GPU_DEVICE AMREX_FORCE_INLINE bool IsSaturated(int i, int j, int k, const amrex::Array4< const amrex::Real > &cell_data, int qc_index)
Definition: ERF_RichardsonNumber.H:80
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real StabilityFunction(amrex::Real Ri, amrex::Real Ri_crit)
Definition: ERF_RichardsonNumber.H:246
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeVerticalShear2(int i, int j, int k, amrex::Real dzInv, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v)
Definition: ERF_RichardsonNumber.H:201
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeRichardson(amrex::Real N2_moist, amrex::Real S2_vert)
Definition: ERF_RichardsonNumber.H:230
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeN2(int i, int j, int k, amrex::Real dzInv, amrex::Real const_grav, const amrex::Array4< const amrex::Real > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_RichardsonNumber.H:126
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ theta
Definition: ERF_MM5.H:20
@ rho
Definition: ERF_Kessler.H:22
@ qc
Definition: ERF_SatAdj.H:36
Definition: ERF_DataStruct.H:99
int qv
Definition: ERF_DataStruct.H:100
int qc
Definition: ERF_DataStruct.H:101