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;
176 const amrex::Array4<const amrex::EBCellFlag>& c_cflag,
179 const amrex::Array4<const amrex::Real>& cell_data,
184 const int qv_index = moisture_indices.
qv;
185 const int qc_index = moisture_indices.
qc;
189 bool saturated =
IsSaturated(i, j, k, cell_data, qc_index);
194 bool is_covered_kp1 = c_cflag(i,j,k+1).isCovered();
195 bool is_covered_km1 = c_cflag(i,j,k-1).isCovered();
197 if (saturated && qc_index >= 0) {
202 inv_theta = 1.0 / theta_l;
204 if (is_covered_kp1) {
206 dtheta_dz = (3.*theta_l -4.*theta_l_km1 + theta_l_km2) * 0.5 * dzInv;
207 }
else if (is_covered_km1) {
209 dtheta_dz = (-theta_l_kp2 + 4.*theta_l_kp1 - 3.*theta_l) * 0.5 * dzInv;
211 dtheta_dz = 0.5 * (theta_l_kp1 - theta_l_km1) * dzInv;
214 }
else if (qv_index >= 0) {
219 inv_theta = 1.0 / theta_v;
221 if (is_covered_kp1) {
223 dtheta_dz = (3.*theta_v -4.*theta_v_km1 + theta_v_km2) * 0.5 * dzInv;
224 }
else if (is_covered_km1) {
226 dtheta_dz = (-theta_v_kp2 + 4.*theta_v_kp1 - 3.*theta_v) * 0.5 * dzInv;
228 dtheta_dz = 0.5 * (theta_v_kp1 - theta_v_km1) * dzInv;
233 amrex::Real theta_kp1 = cell_data(i, j, k+1, rhoTheta_comp) / cell_data(i, j, k+1, rho_comp);
234 amrex::Real theta_km1 = cell_data(i, j, k-1, rhoTheta_comp) / cell_data(i, j, k-1, rho_comp);
236 inv_theta = 1.0 /
theta;
238 if (is_covered_kp1) {
239 amrex::Real theta_km2 = cell_data(i, j, k-2, rhoTheta_comp) / cell_data(i, j, k-2, rho_comp);
240 dtheta_dz = (3.*
theta - 4.*theta_km1 + theta_km2) * 0.5 * dzInv;
241 }
else if (is_covered_km1) {
242 amrex::Real theta_kp2 = cell_data(i, j, k+2, rhoTheta_comp) / cell_data(i, j, k+2, rho_comp);
243 dtheta_dz = (-theta_kp2 + 4.*theta_kp1 - 3.*
theta) * 0.5 * dzInv;
245 dtheta_dz = 0.5 * (theta_kp1 - theta_km1) * dzInv;
249 return const_grav * inv_theta * dtheta_dz;
283 const amrex::Array4<const amrex::Real>& u,
284 const amrex::Array4<const amrex::Real>& v)
286 amrex::Real u_c_km1 = 0.5 * (u(i, j, k-1) + u(i+1, j, k-1));
287 amrex::Real u_c_kp1 = 0.5 * (u(i, j, k+1) + u(i+1, j, k+1));
289 amrex::Real v_c_km1 = 0.5 * (v(i, j, k-1) + v(i, j+1, k-1));
290 amrex::Real v_c_kp1 = 0.5 * (v(i, j, k+1) + v(i, j+1, k+1));
292 amrex::Real du_dz = 0.5 * (u_c_kp1 - u_c_km1) * dzInv;
293 amrex::Real dv_dz = 0.5 * (v_c_kp1 - v_c_km1) * dzInv;
295 amrex::Real S2_vert = du_dz * du_dz + dv_dz * dv_dz;
304 const amrex::Array4<const amrex::EBCellFlag>& c_cflag,
305 const amrex::Array4<const amrex::Real>& u_vfrac,
306 const amrex::Array4<const amrex::Real>& v_vfrac,
308 const amrex::Array4<const amrex::Real>& u,
309 const amrex::Array4<const amrex::Real>& v)
314 amrex::Real u_c = ( u_vfrac(i,j,k) * u(i, j, k) + u_vfrac(i+1,j,k) * u(i+1, j, k) )
315 / ( u_vfrac(i,j,k) + u_vfrac(i+1,j,k) );
316 amrex::Real v_c = ( v_vfrac(i,j,k) * v(i, j, k) + v_vfrac(i,j+1,k) * v(i, j+1, k) )
317 / ( v_vfrac(i,j,k) + v_vfrac(i,j+1,k) );
319 if (c_cflag(i,j,k+1).isCovered()) {
320 amrex::Real u_c_km1 = ( u_vfrac(i,j,k-1) * u(i, j, k-1) + u_vfrac(i+1,j,k-1) * u(i+1, j, k-1) )
321 / ( u_vfrac(i,j,k-1) + u_vfrac(i+1,j,k-1) );
322 amrex::Real u_c_km2 = ( u_vfrac(i,j,k-2) * u(i, j, k-2) + u_vfrac(i+1,j,k-2) * u(i+1, j, k-2) )
323 / ( u_vfrac(i,j,k-2) + u_vfrac(i+1,j,k-2) );
324 amrex::Real v_c_km1 = ( v_vfrac(i,j,k-1) * v(i, j, k-1) + v_vfrac(i,j+1,k-1) * v(i, j+1, k-1) )
325 / ( v_vfrac(i,j,k-1) + v_vfrac(i,j+1,k-1) );
326 amrex::Real v_c_km2 = ( v_vfrac(i,j,k-2) * v(i, j, k-2) + v_vfrac(i,j+1,k-2) * v(i, j+1, k-2) )
327 / ( v_vfrac(i,j,k-2) + v_vfrac(i,j+1,k-2) );
328 du_dz = (3.*u_c - 4.*u_c_km1 + u_c_km2) * 0.5 * dzInv;
329 dv_dz = (3.*v_c - 4.*v_c_km1 + v_c_km2) * 0.5 * dzInv;
331 }
else if (c_cflag(i,j,k-1).isCovered()) {
332 amrex::Real u_c_kp1 = ( u_vfrac(i,j,k+1) * u(i, j, k+1) + u_vfrac(i+1,j,k+1) * u(i+1, j, k+1) )
333 / ( u_vfrac(i,j,k+1) + u_vfrac(i+1,j,k+1) );
334 amrex::Real u_c_kp2 = ( u_vfrac(i,j,k+2) * u(i, j, k+2) + u_vfrac(i+1,j,k+2) * u(i+1, j, k+2) )
335 / ( u_vfrac(i,j,k+2) + u_vfrac(i+1,j,k+2) );
336 amrex::Real v_c_kp1 = ( v_vfrac(i,j,k+1) * v(i, j, k+1) + v_vfrac(i,j+1,k+1) * v(i, j+1, k+1) )
337 / ( v_vfrac(i,j,k+1) + v_vfrac(i,j+1,k+1) );
338 amrex::Real v_c_kp2 = ( v_vfrac(i,j,k+2) * v(i, j, k+2) + v_vfrac(i,j+1,k+2) * v(i, j+1, k+2) )
339 / ( v_vfrac(i,j,k+2) + v_vfrac(i,j+1,k+2) );
340 du_dz = (-u_c_kp2 + 4.*u_c_kp1 - 3.*u_c) * 0.5 * dzInv;
341 dv_dz = (-v_c_kp2 + 4.*v_c_kp1 - 3.*v_c) * 0.5 * dzInv;
343 amrex::Real u_c_km1 = ( u_vfrac(i,j,k-1) * u(i, j, k-1) + u_vfrac(i+1,j,k-1) * u(i+1, j, k-1) )
344 / ( u_vfrac(i,j,k-1) + u_vfrac(i+1,j,k-1) );
345 amrex::Real u_c_kp1 = ( u_vfrac(i,j,k+1) * u(i, j, k+1) + u_vfrac(i+1,j,k+1) * u(i+1, j, k+1) )
346 / ( u_vfrac(i,j,k+1) + u_vfrac(i+1,j,k+1) );
347 amrex::Real v_c_km1 = ( v_vfrac(i,j,k-1) * v(i, j, k-1) + v_vfrac(i,j+1,k-1) * v(i, j+1, k-1) )
348 / ( v_vfrac(i,j,k-1) + v_vfrac(i,j+1,k-1) );
349 amrex::Real v_c_kp1 = ( v_vfrac(i,j,k+1) * v(i, j, k+1) + v_vfrac(i,j+1,k+1) * v(i, j+1, k+1) )
350 / ( v_vfrac(i,j,k+1) + v_vfrac(i,j+1,k+1) );
351 du_dz = 0.5 * (u_c_kp1 - u_c_km1) * dzInv;
352 dv_dz = 0.5 * (v_c_kp1 - v_c_km1) * dzInv;
355 amrex::Real S2_vert = du_dz * du_dz + dv_dz * dv_dz;
372 amrex::Real S2_safe = amrex::max(S2_vert, 1.0e-10);
373 return N2_moist / S2_safe;
390 }
else if (Ri > 0.0) {
391 return std::sqrt(1.0 - Ri / Ri_crit);
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
rho
Definition: ERF_InitCustomPert_Bubble.H:106
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:386
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:281
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeVerticalShear2_EB(int i, int j, int k, const amrex::Array4< const amrex::EBCellFlag > &c_cflag, const amrex::Array4< const amrex::Real > &u_vfrac, const amrex::Array4< const amrex::Real > &v_vfrac, amrex::Real dzInv, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &v)
Definition: ERF_RichardsonNumber.H:303
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeRichardson(amrex::Real N2_moist, amrex::Real S2_vert)
Definition: ERF_RichardsonNumber.H:370
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_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeN2_EB(int i, int j, int k, const amrex::Array4< const amrex::EBCellFlag > &c_cflag, amrex::Real dzInv, amrex::Real const_grav, const amrex::Array4< const amrex::Real > &cell_data, const MoistureComponentIndices &moisture_indices)
Definition: ERF_RichardsonNumber.H:175
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ theta
Definition: ERF_MM5.H:20
@ qc
Definition: ERF_SatAdj.H:36
Definition: ERF_DataStruct.H:106
int qv
Definition: ERF_DataStruct.H:107
int qc
Definition: ERF_DataStruct.H:108