ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_RichardsonNumber.H
Go to the documentation of this file.
1 /** \file ERF_RichardsonNumber.H
2  *
3  * Functions for computing moist Richardson number and Mason (1989) stability function
4  * for Smagorinsky turbulence closure with conditional instability correction.
5  *
6  * PHYSICAL BACKGROUND:
7  * -------------------
8  * In cloud-resolving simulations, standard Smagorinsky models overpredict turbulent
9  * mixing in conditionally unstable regions (cloud updrafts) where buoyancy-driven
10  * instability occurs. This module implements the moist Richardson number correction
11  * following Mason (1989) and Stevens et al. (2005) for marine stratocumulus.
12  *
13  * KEY EQUATIONS:
14  * --------------
15  * 1. Moist Brunt-Väisälä frequency:
16  * N²_moist = (g/θ_v) * ∂θ_v/∂z (unsaturated)
17  * N²_moist = (g/θ_l) * ∂θ_l/∂z (saturated, in clouds)
18  *
19  * 2. Moist Richardson number:
20  * Ri_moist = N²_moist / S²_vert
21  * where S²_vert = (∂u/∂z)² + (∂v/∂z)²
22  *
23  * 3. Mason (1989) stability function:
24  * f(Ri) = 0 if Ri ≥ Ri_crit (complete suppression)
25  * = sqrt(1 - Ri/Ri_crit) if 0 < Ri < Ri_crit (partial suppression)
26  * = 1 if Ri ≤ 0 (no suppression, unstable/neutral)
27  *
28  * 4. Corrected vertical eddy viscosity:
29  * μ_t,v = ρ(C_s·Δ)²|S| · f(Ri_moist)
30  *
31  * PHYSICAL INTERPRETATION:
32  * ------------------------
33  * - In cloud updrafts: ∂θ_l/∂z < 0 → Ri < 0 → f(Ri) = 1.0 → full mixing (correct!)
34  * - In stable regions: ∂θ_v/∂z > 0 → Ri > 0 → f(Ri) < 1.0 → reduced mixing
35  * - In neutral regions: Ri ≈ 0 → f(Ri) ≈ 1.0 → standard Smagorinsky
36  *
37  * CUDA-CONSERVATIVE DESIGN:
38  * --------------------------
39  * - All functions are explicit device functions (not lambdas)
40  * - Use simple scalar parameters (no struct captures in device context)
41  * - Geometry passed via local scalars (dzInv), not complex objects
42  * - Reuse existing ERF GetThetav/GetThetal from ERF_MoistUtils.H
43  *
44  * REFERENCES:
45  * -----------
46  * - Mason, P. J. (1989): Large-eddy simulation of the convective atmospheric boundary
47  * layer. J. Atmos. Sci., 46, 1492-1516.
48  * - Stevens et al. (2005): Evaluation of large-eddy simulations via observations of
49  * nocturnal marine stratocumulus. Mon. Wea. Rev., 133, 1443-1462.
50  * - Smagorinsky, J. (1963): General circulation experiments with the primitive equations.
51  * Mon. Wea. Rev., 91, 99-164.
52  */
53 
54 #ifndef ERF_RICHARDSON_NUMBER_H_
55 #define ERF_RICHARDSON_NUMBER_H_
56 
57 #include <AMReX_GpuQualifiers.H>
58 #include <AMReX_REAL.H>
59 #include <cmath>
60 #include "ERF_Constants.H"
61 #include "ERF_IndexDefines.H"
62 #include "ERF_DataStruct.H"
63 #include "ERF_MoistUtils.H"
64 
65 /**
66  * Check if grid cell is saturated based on cloud liquid water content.
67  *
68  * @param[in] i,j,k cell indices
69  * @param[in] cell_data conserved state (ρ, ρθ, ρq_v, ρq_c, etc.)
70  * @param[in] moisture_indices indices for moisture species
71  * @return true if saturated (q_c > 1e-8 kg/kg), false otherwise
72  *
73  * PHYSICAL JUSTIFICATION: Threshold of 1e-8 kg/kg (0.01 g/kg) is typical
74  * for distinguishing cloudy from clear regions in LES. Below this threshold,
75  * condensate is negligible and air is effectively unsaturated.
76  */
77 AMREX_GPU_DEVICE
78 AMREX_FORCE_INLINE
79 bool
80 IsSaturated(int i, int j, int k,
81  const amrex::Array4<const amrex::Real>& cell_data,
82  int qc_index)
83 {
84  if (qc_index >= 0) {
85  amrex::Real rho = cell_data(i, j, k, Rho_comp);
86  amrex::Real qc = cell_data(i, j, k, qc_index) / rho;
87  return (qc > 1.0e-8); // 0.01 g/kg threshold
88  }
89  return false;
90 }
91 
92 /**
93  * Compute moist Brunt-Väisälä frequency squared (N²_moist).
94  *
95  * @param[in] i,j,k cell indices
96  * @param[in] dzInv inverse vertical grid spacing (1/Δz) [1/m]
97  * @param[in] const_grav gravitational acceleration g [m/s²]
98  * @param[in] cell_data conserved state
99  * @param[in] moisture_indices moisture species indices
100  * @return N²_moist [1/s²]
101  *
102  * PHYSICS:
103  * --------
104  * Uses centered finite differences: ∂/∂z ≈ (f_{k+1} - f_{k-1}) / (2Δz)
105  *
106  * Unsaturated: N² = (g/θ_v) · ∂θ_v/∂z
107  * where θ_v = θ(1 + 0.61q_v - q_c) accounts for vapor buoyancy and condensate loading
108  *
109  * Saturated: N² = (g/θ_l) · ∂θ_l/∂z
110  * where θ_l = θ - (L_v/c_p·Π)·q_c is the liquid water potential temperature
111  *
112  * CONDITIONAL INSTABILITY:
113  * ------------------------
114  * In cloud updrafts, ∂θ_l/∂z can be **negative** due to release of latent heat
115  * during ascent. This gives N² < 0, which is physically correct for convective
116  * instability and should enhance (not suppress) mixing.
117  *
118  * UNITS:
119  * ------
120  * Input: dzInv [1/m], const_grav [m/s²], θ [K]
121  * Output: N² [1/s²]
122  */
123 AMREX_GPU_DEVICE
124 AMREX_FORCE_INLINE
126 ComputeN2(int i, int j, int k,
127  amrex::Real dzInv,
128  amrex::Real const_grav,
129  const amrex::Array4<const amrex::Real>& cell_data,
130  const MoistureComponentIndices& moisture_indices)
131 {
132  const int rho_comp = Rho_comp;
133  const int rhoTheta_comp = RhoTheta_comp;
134  const int qv_index = moisture_indices.qv;
135  const int qc_index = moisture_indices.qc;
136 
137  amrex::Real rho = cell_data(i, j, k, rho_comp);
138 
139  bool saturated = IsSaturated(i, j, k, cell_data, qc_index);
140 
141  amrex::Real inv_theta;
142  amrex::Real dtheta_dz;
143 
144  if (saturated && qc_index >= 0) {
145  amrex::Real theta_l = GetThetal(i, j, k , cell_data, moisture_indices);
146  amrex::Real theta_l_kp1= GetThetal(i, j, k+1, cell_data, moisture_indices);
147  amrex::Real theta_l_km1= GetThetal(i, j, k-1, cell_data, moisture_indices);
148 
149  inv_theta = 1.0 / theta_l;
150  dtheta_dz = 0.5 * (theta_l_kp1 - theta_l_km1) * dzInv;
151 
152  } else if (qv_index >= 0) {
153  amrex::Real theta_v = GetThetav(i, j, k , cell_data, moisture_indices);
154  amrex::Real theta_v_kp1= GetThetav(i, j, k+1, cell_data, moisture_indices);
155  amrex::Real theta_v_km1= GetThetav(i, j, k-1, cell_data, moisture_indices);
156 
157  inv_theta = 1.0 / theta_v;
158  dtheta_dz = 0.5 * (theta_v_kp1 - theta_v_km1) * dzInv;
159 
160  } else {
161  amrex::Real theta = cell_data(i, j, k, rhoTheta_comp) / rho;
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);
164 
165  inv_theta = 1.0 / theta;
166  dtheta_dz = 0.5 * (theta_kp1 - theta_km1) * dzInv;
167  }
168 
169  return const_grav * inv_theta * dtheta_dz;
170 }
171 
172 /**
173  * Compute vertical shear squared (S²_vert).
174  *
175  * @param[in] i,j,k cell indices
176  * @param[in] dzInv inverse vertical grid spacing (1/Δz) [1/m]
177  * @param[in] u face-centered x-velocity (u) array with ngrow≥1 in z
178  * @param[in] v face-centered y-velocity (v) array with ngrow≥1 in z
179  * @return S²_vert = (∂u/∂z)² + (∂v/∂z)² [1/s²]
180  *
181  * DISCRETIZATION:
182  * ---------------
183  * ERF uses Arakawa C-grid:
184  * u defined at (i+1/2, j, k) (x-faces)
185  * v defined at (i, j+1/2, k) (y-faces)
186  * w defined at (i, j, k+1/2) (z-faces)
187  * cell centers at (i, j, k)
188  *
189  * To get shear at cell center (i,j,k), interpolate velocities to center then difference:
190  * ∂u/∂z|_{i,j,k} ≈ (u_{i,j,k+1} - u_{i,j,k-1}) / (2Δz)
191  * where u_{i,j,k} = 0.5 * (u_{i-1/2,j,k} + u_{i+1/2,j,k})
192  *
193  * UNITS:
194  * ------
195  * Input: dzInv [1/m], u, v [m/s]
196  * Output: S²_vert [1/s²]
197  */
198 AMREX_GPU_DEVICE
199 AMREX_FORCE_INLINE
201 ComputeVerticalShear2(int i, int j, int k,
202  amrex::Real dzInv,
203  const amrex::Array4<const amrex::Real>& u,
204  const amrex::Array4<const amrex::Real>& v)
205 {
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));
208 
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));
211 
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;
214 
215  amrex::Real S2_vert = du_dz * du_dz + dv_dz * dv_dz;
216 
217  return S2_vert;
218 }
219 
220 /**
221  * Compute moist Richardson number.
222  *
223  * @param[in] N2_moist moist Brunt-Väisälä frequency squared [1/s²]
224  * @param[in] S2_vert vertical wind shear squared [1/s²]
225  * @return Ri_moist = N²_moist / S²_vert [dimensionless]
226  */
227 AMREX_GPU_DEVICE
228 AMREX_FORCE_INLINE
231 {
232  amrex::Real S2_safe = amrex::max(S2_vert, 1.0e-10);
233  return N2_moist / S2_safe;
234 }
235 
236 /**
237  * Mason (1989) stability function for turbulent mixing suppression.
238  *
239  * @param[in] Ri gradient Richardson number [dimensionless]
240  * @param[in] Ri_crit critical Richardson number (default 0.25)
241  * @return f(Ri) ∈ [0, 1] [dimensionless]
242  */
243 AMREX_GPU_DEVICE
244 AMREX_FORCE_INLINE
247 {
248  if (Ri >= Ri_crit) {
249  return 0.0;
250  } else if (Ri > 0.0) {
251  return std::sqrt(1.0 - Ri / Ri_crit);
252  } else {
253  return 1.0;
254  }
255 }
256 
257 #endif // ERF_RICHARDSON_NUMBER_H_
#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