ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_MakeMomSources.cpp File Reference
#include <AMReX_MultiFab.H>
#include <AMReX_ArrayLim.H>
#include <AMReX_BCRec.H>
#include <AMReX_TableData.H>
#include <AMReX_GpuContainers.H>
#include "ERF_NumericalDiffusion.H"
#include "ERF_PlaneAverage.H"
#include "ERF_TI_slow_headers.H"
#include "ERF_SrcHeaders.H"
#include "ERF_Utils.H"
Include dependency graph for ERF_MakeMomSources.cpp:

Functions

void make_mom_sources (Real time, Real, const Vector< MultiFab > &S_data, const MultiFab *z_phys_nd, const MultiFab *z_phys_cc, Vector< Real > &stretched_dz_h, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &wvel, MultiFab &xmom_src, MultiFab &ymom_src, MultiFab &zmom_src, const MultiFab &base_state, MultiFab *forest_drag, MultiFab *terrain_blank, MultiFab *cosPhi_mf, MultiFab *sinPhi_mf, const Geometry geom, const SolverChoice &solverChoice, Vector< std::unique_ptr< MultiFab >> &, const Real *dptr_u_geos, const Real *dptr_v_geos, const Real *dptr_wbar_sub, const Vector< Real * > d_rayleigh_ptrs_at_lev, const amrex::Real *d_sinesq_at_lev, const amrex::Real *d_sinesq_stag_at_lev, const Vector< Real * > d_sponge_ptrs_at_lev, const Vector< MultiFab > *forecast_state_at_lev, const MultiFab *surface_state_at_lev, InputSoundingData &input_sounding_data, const eb_ &ebfact, bool is_slow_step)
 

Function Documentation

◆ make_mom_sources()

void make_mom_sources ( Real  time,
Real  ,
const Vector< MultiFab > &  S_data,
const MultiFab *  z_phys_nd,
const MultiFab *  z_phys_cc,
Vector< Real > &  stretched_dz_h,
const MultiFab &  xvel,
const MultiFab &  yvel,
const MultiFab &  wvel,
MultiFab &  xmom_src,
MultiFab &  ymom_src,
MultiFab &  zmom_src,
const MultiFab &  base_state,
MultiFab *  forest_drag,
MultiFab *  terrain_blank,
MultiFab *  cosPhi_mf,
MultiFab *  sinPhi_mf,
const Geometry  geom,
const SolverChoice solverChoice,
Vector< std::unique_ptr< MultiFab >> &  ,
const Real dptr_u_geos,
const Real dptr_v_geos,
const Real dptr_wbar_sub,
const Vector< Real * >  d_rayleigh_ptrs_at_lev,
const amrex::Real d_sinesq_at_lev,
const amrex::Real d_sinesq_stag_at_lev,
const Vector< Real * >  d_sponge_ptrs_at_lev,
const Vector< MultiFab > *  forecast_state_at_lev,
const MultiFab *  surface_state_at_lev,
InputSoundingData input_sounding_data,
const eb_ ebfact,
bool  is_slow_step 
)

Function for computing the slow RHS for the evolution equations for the density, potential temperature and momentum.

Parameters
[in]timecurrent time
[in]dtcurrent slow or fast timestep size
[in]S_datacurrent solution
[in]xvelx-component of velocity
[in]yvely-component of velocity
[in]xmom_srcsource terms for x-momentum
[in]ymom_srcsource terms for y-momentum
[in]zmom_srcsource terms for z-momentum
[in]geomContainer for geometric information
[in]solverChoiceContainer for solver parameters
[in]mapfacmap factors
[in]dptr_u_geoscustom geostrophic wind profile
[in]dptr_v_geoscustom geostrophic wind profile
[in]dptr_wbar_subsubsidence source term
[in]d_rayleigh_ptrs_at_levVector of {strength of Rayleigh damping, reference value for xvel/yvel/zvel/theta} used to define Rayleigh damping
[in]d_sinesq_at_levsin( (pi/2) (z-z_t)/(damping depth)) at cell centers
[in]d_sinesq_stag_at_levsin( (pi/2) (z-z_t)/(damping depth)) at z-faces
69 {
70  BL_PROFILE_REGION("erf_make_mom_sources()");
71 
72  Box domain(geom.Domain());
73  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
74 
75  // Initialize sources to zero each time we may use them
76  xmom_src.setVal(0.0);
77  ymom_src.setVal(0.0);
78  zmom_src.setVal(0.0);
79 
80  MultiFab r_hse (base_state, make_alias, BaseState::r0_comp , 1);
81 
82  // flags to apply certain source terms in substep call only
83  bool use_Rayleigh_fast_uv = ( (solverChoice.dampingChoice.rayleigh_damping_type == RayleighDampingType::FastExplicit) ||
85  bool use_Rayleigh_fast_w = (solverChoice.dampingChoice.rayleigh_damping_type == RayleighDampingType::FastExplicit);
86  bool use_canopy_fast = solverChoice.forest_substep;
87  bool use_ImmersedForcing_fast = solverChoice.immersed_forcing_substep;
88 
89  // *****************************************************************************
90  // Define source term for all three components of momenta from
91  // one Coriolis forcing for (xmom,ymom,zmom)
92  // two Rayleigh damping for (xmom,ymom,zmom)
93  // three Constant / height-dependent geostrophic forcing
94  // Real(4.) Subsidence
95  // Real(5.) Nudging towards input sounding data
96  // Real(6.) Numerical diffusion for (xmom,ymom,zmom)
97  // Real(7.) Sponge
98  // Real(8.) Forest canopy
99  // 9a. Immersed forcing for terrain
100  // 9b. Immersed forcing for buildings
101  // Real(10.) Constant mass flux
102  // *****************************************************************************
103  // NOTE: buoyancy is now computed in a separate routine - it should not appear here
104  // *****************************************************************************
105  //const bool l_use_ndiff = solverChoice.use_num_diff;
106 
107  if (solverChoice.terrain_type == TerrainType::ImmersedForcing) {
108  if (solverChoice.do_forest_drag) {
109  amrex::Error(" Currently forest canopy cannot be used with immersed forcing");
110  }
111  }
112 
113 
114  // *****************************************************************************
115  // Data for Coriolis forcing
116  // *****************************************************************************
117  auto use_coriolis = solverChoice.use_coriolis;
118  auto coriolis_factor = solverChoice.coriolis_factor;
119  auto cosphi = solverChoice.cosphi;
120  auto sinphi = solverChoice.sinphi;
121  auto var_coriolis = solverChoice.variable_coriolis;
122 
123  // *****************************************************************************
124  // Flag for Geostrophic forcing
125  // *****************************************************************************
126  auto abl_geo_forcing = solverChoice.abl_geo_forcing;
127  auto geo_wind_profile = solverChoice.have_geo_wind_profile;
128 
129  // *****************************************************************************
130  // Data for Rayleigh damping
131  // *****************************************************************************
132  auto rayleigh_damp_U = solverChoice.dampingChoice.rayleigh_damp_U;
133  auto rayleigh_damp_V = solverChoice.dampingChoice.rayleigh_damp_V;
134  auto rayleigh_damp_W = solverChoice.dampingChoice.rayleigh_damp_W;
135 
136  Real* ubar = d_rayleigh_ptrs_at_lev[Rayleigh::ubar];
137  Real* vbar = d_rayleigh_ptrs_at_lev[Rayleigh::vbar];
138  Real* wbar = d_rayleigh_ptrs_at_lev[Rayleigh::wbar];
139 
140  // *****************************************************************************
141  // Data for constant mass flux
142  // *****************************************************************************
143  bool enforce_massflux_x = (solverChoice.const_massflux_u != 0);
144  bool enforce_massflux_y = (solverChoice.const_massflux_v != 0);
145  Real U_target = solverChoice.const_massflux_u;
146  Real V_target = solverChoice.const_massflux_v;
147  int massflux_klo = solverChoice.massflux_klo;
148  int massflux_khi = solverChoice.massflux_khi;
149 
150  // These will be updated by integrating through the planar average profiles
151  Real rhoUA_target{0};
152  Real rhoVA_target{0};
153  Real rhoUA{0};
154  Real rhoVA{0};
155 
156  // *****************************************************************************
157  // Planar averages for subsidence, nudging, or constant mass flux
158  // *****************************************************************************
159  Table1D<Real> dptr_r_plane, dptr_u_plane, dptr_v_plane;
160  TableData<Real, 1> r_plane_tab, u_plane_tab, v_plane_tab;
161 
162  if (is_slow_step && (dptr_wbar_sub || solverChoice.nudging_from_input_sounding ||
163  enforce_massflux_x || enforce_massflux_y))
164  {
165  // The plane averaging operates at fixed z not fixed height so is not correct for variable dz
166  AMREX_ALWAYS_ASSERT(solverChoice.mesh_type != MeshType::VariableDz);
167 
168  const int offset = 1;
169  const int u_offset = 1;
170  const int v_offset = 1;
171 
172  //
173  // We use the alias here to control ncomp inside the PlaneAverage
174  //
175  MultiFab cons(S_data[IntVars::cons], make_alias, 0, 1);
176 
177  IntVect ng_c = S_data[IntVars::cons].nGrowVect(); ng_c[2] = offset;
178  PlaneAverage r_ave(&cons, geom, solverChoice.ave_plane, ng_c);
179  r_ave.compute_averages(ZDir(), r_ave.field());
180 
181  int ncell = r_ave.ncell_line();
182  Gpu::HostVector< Real> r_plane_h(ncell);
183  Gpu::DeviceVector< Real> r_plane_d(ncell);
184 
185  r_ave.line_average(Rho_comp, r_plane_h);
186 
187  Gpu::copyAsync(Gpu::hostToDevice, r_plane_h.begin(), r_plane_h.end(), r_plane_d.begin());
188 
189  Real* dptr_r = r_plane_d.data();
190 
191  Box tdomain = domain; tdomain.grow(2,ng_c[2]);
192  r_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
193 
194  dptr_r_plane = r_plane_tab.table();
195  ParallelFor(ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
196  {
197  dptr_r_plane(k-offset) = dptr_r[k];
198  });
199 
200  // U and V momentum
201  IntVect ng_u = S_data[IntVars::xmom].nGrowVect(); ng_u[2] = u_offset;
202  PlaneAverage u_ave(&(S_data[IntVars::xmom]), geom, solverChoice.ave_plane, ng_u);
203 
204  IntVect ng_v = S_data[IntVars::ymom].nGrowVect(); ng_v[2] = v_offset;
205  PlaneAverage v_ave(&(S_data[IntVars::ymom]), geom, solverChoice.ave_plane, ng_v);
206 
207  u_ave.compute_averages(ZDir(), u_ave.field());
208  v_ave.compute_averages(ZDir(), v_ave.field());
209 
210  int u_ncell = u_ave.ncell_line();
211  int v_ncell = v_ave.ncell_line();
212  Gpu::HostVector< Real> u_plane_h(u_ncell), v_plane_h(v_ncell);
213  Gpu::DeviceVector< Real> u_plane_d(u_ncell), v_plane_d(v_ncell);
214 
215  u_ave.line_average(0, u_plane_h);
216  v_ave.line_average(0, v_plane_h);
217 
218  Gpu::copyAsync(Gpu::hostToDevice, u_plane_h.begin(), u_plane_h.end(), u_plane_d.begin());
219  Gpu::copyAsync(Gpu::hostToDevice, v_plane_h.begin(), v_plane_h.end(), v_plane_d.begin());
220 
221  Real* dptr_u = u_plane_d.data();
222  Real* dptr_v = v_plane_d.data();
223 
224  Box udomain = domain; udomain.grow(2,ng_u[2]);
225  Box vdomain = domain; vdomain.grow(2,ng_v[2]);
226  u_plane_tab.resize({udomain.smallEnd(2)}, {udomain.bigEnd(2)});
227  v_plane_tab.resize({vdomain.smallEnd(2)}, {vdomain.bigEnd(2)});
228 
229  dptr_u_plane = u_plane_tab.table();
230  ParallelFor(u_ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
231  {
232  dptr_u_plane(k-u_offset) = dptr_u[k];
233  });
234 
235  dptr_v_plane = v_plane_tab.table();
236  ParallelFor(v_ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
237  {
238  dptr_v_plane(k-v_offset) = dptr_v[k];
239  });
240 
241  // sum in z for massflux adjustment
242  if (enforce_massflux_x || enforce_massflux_y) {
243  Real Lx = geom.ProbHi(0) - geom.ProbLo(0);
244  Real Ly = geom.ProbHi(1) - geom.ProbLo(1);
245 
246  if (solverChoice.mesh_type == MeshType::ConstantDz) {
247  // note: massflux_khi corresponds to unstaggered indices in this case
248  rhoUA = std::accumulate(u_plane_h.begin() + u_offset + massflux_klo,
249  u_plane_h.begin() + u_offset + massflux_khi+1, zero);
250  rhoVA = std::accumulate(v_plane_h.begin() + v_offset + massflux_klo,
251  v_plane_h.begin() + v_offset + massflux_khi+1, zero);
252  rhoUA_target = std::accumulate(r_plane_h.begin() + offset + massflux_klo,
253  r_plane_h.begin() + offset + massflux_khi+1, zero);
254  rhoVA_target = rhoUA_target;
255 
256  rhoUA *= geom.CellSize(2) * Ly;
257  rhoVA *= geom.CellSize(2) * Lx;
258  rhoUA_target *= geom.CellSize(2) * Ly;
259  rhoVA_target *= geom.CellSize(2) * Lx;
260 
261  } else if (solverChoice.mesh_type == MeshType::StretchedDz) {
262  // note: massflux_khi corresponds to staggered indices in this case
263  for (int k=massflux_klo; k < massflux_khi; ++k) {
264  rhoUA += u_plane_h[k + u_offset] * stretched_dz_h[k];
265  rhoVA += v_plane_h[k + v_offset] * stretched_dz_h[k];
266  rhoUA_target += r_plane_h[k + offset] * stretched_dz_h[k];
267  }
268  rhoVA_target = rhoUA_target;
269 
270  rhoUA *= Ly;
271  rhoVA *= Lx;
272  rhoUA_target *= Ly;
273  rhoVA_target *= Lx;
274  }
275 
276  // at this point, this is integrated rho*dA
277  rhoUA_target *= U_target;
278  rhoVA_target *= V_target;
279 
280  Print() << "Integrated mass flux : " << rhoUA << " " << rhoVA
281  << " (target: " << rhoUA_target << " " << rhoVA_target << ")"
282  << std::endl;
283  }
284  }
285 
286  // *****************************************************************************
287  // Add all the other forcings
288  // *****************************************************************************
289  for ( MFIter mfi(S_data[IntVars::cons]); mfi.isValid(); ++mfi)
290  {
291  Box tbx = mfi.nodaltilebox(0);
292  Box tby = mfi.nodaltilebox(1);
293  Box tbz = mfi.nodaltilebox(2);
294  if (tbz.bigEnd(2) == domain.bigEnd(2)+1) tbz.growHi(2,-1);
295 
296  const Array4<const Real>& cell_data = S_data[IntVars::cons].array(mfi);
297  const Array4<const Real>& rho_u = S_data[IntVars::xmom].array(mfi);
298  const Array4<const Real>& rho_v = S_data[IntVars::ymom].array(mfi);
299  const Array4<const Real>& rho_w = S_data[IntVars::zmom].array(mfi);
300 
301  const Array4<const Real>& u = xvel.array(mfi);
302  const Array4<const Real>& v = yvel.array(mfi);
303  const Array4<const Real>& w = wvel.array(mfi);
304 
305  const Array4< Real>& xmom_src_arr = xmom_src.array(mfi);
306  const Array4< Real>& ymom_src_arr = ymom_src.array(mfi);
307  const Array4< Real>& zmom_src_arr = zmom_src.array(mfi);
308 
309  const Array4<const Real>& r0 = r_hse.const_array(mfi);
310 
311  const Array4<const Real>& f_drag_arr = (forest_drag) ? forest_drag->const_array(mfi) :
312  Array4<const Real>{};
313  const Array4<const Real>& t_blank_arr = (terrain_blank) ? terrain_blank->const_array(mfi) :
314  Array4<const Real>{};
315 
316  const Array4<const Real>& cphi_arr = (cosPhi_mf) ? cosPhi_mf->const_array(mfi) :
317  Array4<const Real>{};
318  const Array4<const Real>& sphi_arr = (sinPhi_mf) ? sinPhi_mf->const_array(mfi) :
319  Array4<const Real>{};
320 
321  const Array4<const Real>& z_nd_arr = z_phys_nd->const_array(mfi);
322  const Array4<const Real>& z_cc_arr = z_phys_cc->const_array(mfi);
323 
324 
325  // *****************************************************************************
326  // one Add CORIOLIS forcing (this assumes east is +x, north is +y)
327  // *****************************************************************************
328  if (use_coriolis && is_slow_step) {
329  if(solverChoice.init_type == InitType::HindCast) {
330  const Array4<const Real>& latlon_arr = (*forecast_state_at_lev)[4].array(mfi);
331  ParallelFor(tbx, tby, tbz,
332  [=] AMREX_GPU_DEVICE (int i, int j, int k)
333  {
334  Real rho_v_loc = fourth * (rho_v(i,j+1,k) + rho_v(i,j,k) + rho_v(i-1,j+1,k) + rho_v(i-1,j,k));
335  Real rho_w_loc = fourth * (rho_w(i,j,k+1) + rho_w(i,j,k) + rho_w(i,j-1,k+1) + rho_w(i,j-1,k));
336  Real latitude = latlon_arr(i,j,k,0);
337  Real sphi_loc = std::sin(latitude*PI/Real(180.0));
338  Real cphi_loc = std::cos(latitude*PI/Real(180.0));
339  xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sphi_loc - rho_w_loc * cphi_loc);
340  },
341  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
342  Real rho_u_loc = fourth * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j-1,k) + rho_u(i,j-1,k));
343  Real latitude = latlon_arr(i,j,k,0);
344  Real sphi_loc = std::sin(latitude*PI/Real(180.0));
345  ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sphi_loc;
346  },
347  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
348  Real rho_u_loc = fourth * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j,k-1) + rho_u(i,j,k-1));
349  Real latitude = latlon_arr(i,j,k,0);
350  Real cphi_loc = std::cos(latitude*PI/Real(180.0));
351  zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cphi_loc;
352  });
353  }
354  else if (var_coriolis && (sinPhi_mf) && (cosPhi_mf)) {
355  ParallelFor(tbx, tby, tbz,
356  [=] AMREX_GPU_DEVICE (int i, int j, int k)
357  {
358  Real rho_v_loc = fourth * (rho_v(i,j+1,k) + rho_v(i,j,k) + rho_v(i-1,j+1,k) + rho_v(i-1,j,k));
359  Real rho_w_loc = fourth * (rho_w(i,j,k+1) + rho_w(i,j,k) + rho_w(i-1,j,k+1) + rho_w(i-1,j,k));
360  Real sphi_loc = myhalf * (sphi_arr(i,j,0) + sphi_arr(i-1,j,0));
361  Real cphi_loc = myhalf * (cphi_arr(i,j,0) + cphi_arr(i-1,j,0));
362  xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sphi_loc - rho_w_loc * cphi_loc);
363  },
364  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
365  Real rho_u_loc = fourth * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j-1,k) + rho_u(i,j-1,k));
366  Real sphi_loc = myhalf * (sphi_arr(i,j,0) + sphi_arr(i,j-1,0));
367  ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sphi_loc;
368  },
369  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
370  Real rho_u_loc = fourth * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j,k-1) + rho_u(i,j,k-1));
371  zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cphi_arr(i,j,0);
372  });
373  } else {
374  if (solverChoice.terrain_type == TerrainType::EB) {
375  Array4<const Real> u_volfrac = (ebfact.get_u_const_factory())->getVolFrac().const_array(mfi);
376  Array4<const Real> v_volfrac = (ebfact.get_v_const_factory())->getVolFrac().const_array(mfi);
377  Array4<const Real> w_volfrac = (ebfact.get_w_const_factory())->getVolFrac().const_array(mfi);
378  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
379  Real rho_v_loc = 0.0;
380  Real rho_w_loc = 0.0;
381  Real v_vol = v_volfrac(i,j+1,k) + v_volfrac(i,j,k) + v_volfrac(i-1,j+1,k) + v_volfrac(i-1,j,k);
382  Real w_vol = w_volfrac(i,j,k+1) + w_volfrac(i,j,k) + w_volfrac(i-1,j,k+1) + w_volfrac(i-1,j,k);
383  if (v_vol > 0.0) {
384  rho_v_loc = ( v_volfrac(i,j+1,k) * rho_v(i,j+1,k) + v_volfrac(i,j,k) * rho_v(i,j,k)
385  + v_volfrac(i-1,j+1,k) * rho_v(i-1,j+1,k) + v_volfrac(i-1,j,k) * rho_v(i-1,j,k)) / v_vol;
386  }
387  if (w_vol > 0.0) {
388  rho_w_loc = ( w_volfrac(i,j,k+1) * rho_w(i,j,k+1) + w_volfrac(i,j,k) * rho_w(i,j,k)
389  + w_volfrac(i-1,j,k+1) * rho_w(i-1,j,k+1) + w_volfrac(i-1,j,k) * rho_w(i-1,j,k)) / w_vol;
390  }
391  xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sinphi - rho_w_loc * cosphi);
392  });
393  ParallelFor(tby, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
394  Real rho_u_loc = 0.0;
395  Real u_vol = u_volfrac(i+1,j,k) + u_volfrac(i,j,k) + u_volfrac(i+1,j-1,k) + u_volfrac(i,j-1,k);
396  if (u_vol > 0.0) {
397  rho_u_loc = ( u_volfrac(i+1,j,k) * rho_u(i+1,j,k) + u_volfrac(i,j,k) * rho_u(i,j,k)
398  + u_volfrac(i+1,j-1,k) * rho_u(i+1,j-1,k) + u_volfrac(i,j-1,k) * rho_u(i,j-1,k)) / u_vol;
399  }
400  ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sinphi;
401  });
402  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
403  Real rho_u_loc = 0.0;
404  Real u_vol = u_volfrac(i+1,j,k) + u_volfrac(i,j,k) + u_volfrac(i+1,j,k-1) + u_volfrac(i,j,k-1);
405  if (u_vol > 0.0) {
406  rho_u_loc = ( u_volfrac(i+1,j,k) * rho_u(i+1,j,k) + u_volfrac(i,j,k) * rho_u(i,j,k)
407  + u_volfrac(i+1,j,k-1) * rho_u(i+1,j,k-1) + u_volfrac(i,j,k-1) * rho_u(i,j,k-1)) / u_vol;
408  }
409  zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cosphi;
410  });
411  } else {
412  ParallelFor(tbx, tby, tbz,
413  [=] AMREX_GPU_DEVICE (int i, int j, int k)
414  {
415  Real rho_v_loc = fourth * (rho_v(i,j+1,k) + rho_v(i,j,k) + rho_v(i-1,j+1,k) + rho_v(i-1,j,k));
416  Real rho_w_loc = fourth * (rho_w(i,j,k+1) + rho_w(i,j,k) + rho_w(i-1,j,k+1) + rho_w(i-1,j,k));
417  xmom_src_arr(i, j, k) += coriolis_factor * (rho_v_loc * sinphi - rho_w_loc * cosphi);
418  },
419  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
420  Real rho_u_loc = fourth * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j-1,k) + rho_u(i,j-1,k));
421  ymom_src_arr(i, j, k) += -coriolis_factor * rho_u_loc * sinphi;
422  },
423  [=] AMREX_GPU_DEVICE (int i, int j, int k) {
424  Real rho_u_loc = fourth * (rho_u(i+1,j,k) + rho_u(i,j,k) + rho_u(i+1,j,k-1) + rho_u(i,j,k-1));
425  zmom_src_arr(i, j, k) += coriolis_factor * rho_u_loc * cosphi;
426  });
427  }
428  } // var_coriolis
429  } // use_coriolis
430 
431  // *****************************************************************************
432  // two Add RAYLEIGH damping
433  // *****************************************************************************
434  Real dampcoef = solverChoice.dampingChoice.rayleigh_dampcoef;
435 
436  if ( (is_slow_step && !use_Rayleigh_fast_uv) || (!is_slow_step && use_Rayleigh_fast_uv)) {
437  if (rayleigh_damp_U) {
438  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
439  {
440  Real rho_on_u_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i-1,j,k,Rho_comp) );
441  Real uu = rho_u(i,j,k) / rho_on_u_face;
442  Real sinesq = d_sinesq_at_lev[k];
443  xmom_src_arr(i, j, k) -= dampcoef*sinesq * (uu - ubar[k]) * rho_on_u_face;
444  });
445  }
446 
447  if (rayleigh_damp_V) {
448  ParallelFor(tby, [=] AMREX_GPU_DEVICE (int i, int j, int k)
449  {
450  Real rho_on_v_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j-1,k,Rho_comp) );
451  Real vv = rho_v(i,j,k) / rho_on_v_face;
452  Real sinesq = d_sinesq_at_lev[k];
453  ymom_src_arr(i, j, k) -= dampcoef*sinesq * (vv - vbar[k]) * rho_on_v_face;
454  });
455  }
456  } // fast or slow step
457 
458  if ( (is_slow_step && !use_Rayleigh_fast_w) || (!is_slow_step && use_Rayleigh_fast_w)) {
459  if (rayleigh_damp_W) {
460  ParallelFor(tbz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
461  {
462  Real rho_on_w_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp) );
463  Real ww = rho_w(i,j,k) / rho_on_w_face;
464  Real sinesq = d_sinesq_stag_at_lev[k];
465  zmom_src_arr(i, j, k) -= dampcoef*sinesq * (ww - wbar[k]) * rho_on_w_face;
466  });
467  }
468  } // fast or slow step
469 
470  // *****************************************************************************
471  // 3a. Add constant GEOSTROPHIC forcing
472  // *****************************************************************************
473  if (is_slow_step) {
474  ParallelFor(tbx, tby, tbz,
475  [=] AMREX_GPU_DEVICE (int i, int j, int k)
476  {
477  Real rho_on_u_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i-1,j,k,Rho_comp) );
478  xmom_src_arr(i, j, k) += rho_on_u_face * abl_geo_forcing[0];
479  },
480  [=] AMREX_GPU_DEVICE (int i, int j, int k)
481  {
482  Real rho_on_v_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j-1,k,Rho_comp) );
483  ymom_src_arr(i, j, k) += rho_on_v_face * abl_geo_forcing[1];
484  },
485  [=] AMREX_GPU_DEVICE (int i, int j, int k)
486  {
487  Real rho_on_w_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp) );
488  zmom_src_arr(i, j, k) += rho_on_w_face * abl_geo_forcing[2];
489  });
490  }
491 
492  // *****************************************************************************
493  // 3b. Add height-dependent GEOSTROPHIC forcing
494  // *****************************************************************************
495  if (geo_wind_profile && is_slow_step) {
496  ParallelFor(tbx, tby,
497  [=] AMREX_GPU_DEVICE (int i, int j, int k)
498  {
499  Real rho_on_u_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i-1,j,k,Rho_comp) );
500  xmom_src_arr(i, j, k) -= coriolis_factor * rho_on_u_face * dptr_v_geos[k] * sinphi;
501  },
502  [=] AMREX_GPU_DEVICE (int i, int j, int k)
503  {
504  Real rho_on_v_face = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j-1,k,Rho_comp) );
505  ymom_src_arr(i, j, k) += coriolis_factor * rho_on_v_face * dptr_u_geos[k] * sinphi;
506  });
507  } // geo_wind_profile
508 
509  // *****************************************************************************
510  // Real(4.) Add custom SUBSIDENCE terms
511  // *****************************************************************************
512  if (solverChoice.custom_w_subsidence && is_slow_step && solverChoice.do_mom_advection) {
513  if (solverChoice.custom_forcing_prim_vars) {
514  const int nr = Rho_comp;
515  ParallelFor(tbx, tby,
516  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
517  {
518  Real dzInv = myhalf*dxInv[2];
519  if (z_nd_arr) {
520  Real z_xf_lo = fourth * ( z_nd_arr(i,j,k ) + z_nd_arr(i,j+1,k )
521  + z_nd_arr(i,j,k-1) + z_nd_arr(i,j+1,k-1) );
522  Real z_xf_hi = fourth * ( z_nd_arr(i,j,k+1) + z_nd_arr(i,j+1,k+1)
523  + z_nd_arr(i,j,k+2) + z_nd_arr(i,j+1,k+2) );
524  dzInv = one / (z_xf_hi - z_xf_lo);
525  }
526  Real rho_on_u_face = myhalf * ( cell_data(i,j,k,nr) + cell_data(i-1,j,k,nr) );
527  Real U_hi = dptr_u_plane(k+1) / dptr_r_plane(k+1);
528  Real U_lo = dptr_u_plane(k-1) / dptr_r_plane(k-1);
529  Real wbar_xf = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
530  xmom_src_arr(i, j, k) -= rho_on_u_face * wbar_xf * (U_hi - U_lo) * dzInv;
531  },
532  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
533  {
534  Real dzInv = myhalf*dxInv[2];
535  if (z_nd_arr) {
536  Real z_yf_lo = fourth * ( z_nd_arr(i,j,k ) + z_nd_arr(i+1,j,k )
537  + z_nd_arr(i,j,k-1) + z_nd_arr(i+1,j,k-1) );
538  Real z_yf_hi = fourth * ( z_nd_arr(i,j,k+1) + z_nd_arr(i+1,j,k+1)
539  + z_nd_arr(i,j,k+2) + z_nd_arr(i+1,j,k+2) );
540  dzInv = one / (z_yf_hi - z_yf_lo);
541  }
542  Real rho_on_v_face = myhalf * ( cell_data(i,j,k,nr) + cell_data(i,j-1,k,nr) );
543  Real V_hi = dptr_v_plane(k+1) / dptr_r_plane(k+1);
544  Real V_lo = dptr_v_plane(k-1) / dptr_r_plane(k-1);
545  Real wbar_yf = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
546  ymom_src_arr(i, j, k) -= rho_on_v_face * wbar_yf * (V_hi - V_lo) * dzInv;
547  });
548  } else {
549  ParallelFor(tbx, tby,
550  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
551  {
552  Real dzInv = myhalf*dxInv[2];
553  if (z_nd_arr) {
554  Real z_xf_lo = fourth * ( z_nd_arr(i,j,k ) + z_nd_arr(i,j+1,k )
555  + z_nd_arr(i,j,k-1) + z_nd_arr(i,j+1,k-1) );
556  Real z_xf_hi = fourth * ( z_nd_arr(i,j,k+1) + z_nd_arr(i,j+1,k+1)
557  + z_nd_arr(i,j,k+2) + z_nd_arr(i,j+1,k+2) );
558  dzInv = one / (z_xf_hi - z_xf_lo);
559  }
560  Real U_hi = dptr_u_plane(k+1) / dptr_r_plane(k+1);
561  Real U_lo = dptr_u_plane(k-1) / dptr_r_plane(k-1);
562  Real wbar_xf = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
563  xmom_src_arr(i, j, k) -= wbar_xf * (U_hi - U_lo) * dzInv;
564  },
565  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
566  {
567  Real dzInv = myhalf*dxInv[2];
568  if (z_nd_arr) {
569  Real z_yf_lo = fourth * ( z_nd_arr(i,j,k ) + z_nd_arr(i+1,j,k )
570  + z_nd_arr(i,j,k-1) + z_nd_arr(i+1,j,k-1) );
571  Real z_yf_hi = fourth * ( z_nd_arr(i,j,k+1) + z_nd_arr(i+1,j,k+1)
572  + z_nd_arr(i,j,k+2) + z_nd_arr(i+1,j,k+2) );
573  dzInv = one / (z_yf_hi - z_yf_lo);
574  }
575  Real V_hi = dptr_v_plane(k+1) / dptr_r_plane(k+1);
576  Real V_lo = dptr_v_plane(k-1) / dptr_r_plane(k-1);
577  Real wbar_yf = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
578  ymom_src_arr(i, j, k) -= wbar_yf * (V_hi - V_lo) * dzInv;
579  });
580  }
581  }
582 
583  // *************************************************************************************
584  // Real(5.) Add nudging towards value specified in input sounding
585  // *************************************************************************************
586  if (solverChoice.nudging_from_input_sounding && is_slow_step)
587  {
588  int itime_n = 0;
589  int itime_np1 = 0;
590  Real coeff_n = one;
591  Real coeff_np1 = zero;
592 
593  Real tau_inv = one / input_sounding_data.tau_nudging;
594 
595  int n_sounding_times = input_sounding_data.input_sounding_time.size();
596 
597  for (int nt = 1; nt < n_sounding_times; nt++) {
598  if (time > input_sounding_data.input_sounding_time[nt]) itime_n = nt;
599  }
600  if (itime_n == n_sounding_times-1) {
601  itime_np1 = itime_n;
602  } else {
603  itime_np1 = itime_n+1;
604  coeff_np1 = (time - input_sounding_data.input_sounding_time[itime_n]) /
605  (input_sounding_data.input_sounding_time[itime_np1] - input_sounding_data.input_sounding_time[itime_n]);
606  coeff_n = one - coeff_np1;
607  }
608 
609  int nr = Rho_comp;
610 
611  const Real* u_inp_sound_n = input_sounding_data.U_inp_sound_d[itime_n].dataPtr();
612  const Real* u_inp_sound_np1 = input_sounding_data.U_inp_sound_d[itime_np1].dataPtr();
613  const Real* v_inp_sound_n = input_sounding_data.V_inp_sound_d[itime_n].dataPtr();
614  const Real* v_inp_sound_np1 = input_sounding_data.V_inp_sound_d[itime_np1].dataPtr();
615  ParallelFor(tbx, tby,
616  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
617  {
618  Real nudge_u = (coeff_n*u_inp_sound_n[k] + coeff_np1*u_inp_sound_np1[k]) - (dptr_u_plane(k)/dptr_r_plane(k));
619  nudge_u *= tau_inv;
620  xmom_src_arr(i, j, k) += cell_data(i, j, k, nr) * nudge_u;
621  },
622  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
623  {
624  Real nudge_v = (coeff_n*v_inp_sound_n[k] + coeff_np1*v_inp_sound_np1[k]) - (dptr_v_plane(k)/dptr_r_plane(k));
625  nudge_v *= tau_inv;
626  ymom_src_arr(i, j, k) += cell_data(i, j, k, nr) * nudge_v;
627  });
628  }
629 
630  // *****************************************************************************
631  // Real(6.) Add NUMERICAL DIFFUSION terms
632  // *****************************************************************************
633 #if 0
634  if (l_use_ndiff) {
635  const Array4<const Real>& mf_ux = mapfac[MapFac::ux]->const_array(mfi);
636  const Array4<const Real>& mf_uy = mapfac[MapFac::uy]->const_array(mfi);
637  const Array4<const Real>& mf_vx = mapfac[MapFac::vx]->const_array(mfi);
638  const Array4<const Real>& mf_vy = mapfac[MapFac::vy]->const_array(mfi);
639  NumericalDiffusion_Xmom(tbx, dt, solverChoice.num_diff_coeff,
640  u, cell_data, xmom_src_arr, mf_ux, mf_uy);
641  NumericalDiffusion_Ymom(tby, dt, solverChoice.num_diff_coeff,
642  v, cell_data, ymom_src_arr, mf_vx, mf_vy);
643  }
644 #endif
645 
646  // *****************************************************************************
647  // Real(7.) Add SPONGING
648  // *****************************************************************************
649  if (is_slow_step) {
650  if (solverChoice.spongeChoice.sponge_type == SpongeType::Input_Sponge)
651  {
652  ApplySpongeZoneBCsForMom_ReadFromFile(solverChoice.spongeChoice, geom, tbx, tby, cell_data,
653  z_cc_arr, xmom_src_arr, ymom_src_arr,
654  rho_u, rho_v, d_sponge_ptrs_at_lev);
655  }
656  else
657  {
658  ApplySpongeZoneBCsForMom(solverChoice.spongeChoice, geom, tbx, tby, tbz,
659  xmom_src_arr, ymom_src_arr, zmom_src_arr, rho_u, rho_v, rho_w,
660  r0, z_nd_arr, z_cc_arr);
661  }
662 
663  if(solverChoice.init_type == InitType::HindCast and solverChoice.hindcast_lateral_forcing){
664 
665  const Array4<const Real>& rho_u_forecast_state = (*forecast_state_at_lev)[IntVars::xmom].array(mfi);
666  const Array4<const Real>& rho_v_forecast_state = (*forecast_state_at_lev)[IntVars::ymom].array(mfi);
667  const Array4<const Real>& rho_w_forecast_state = (*forecast_state_at_lev)[IntVars::zmom].array(mfi);
668  const Array4<const Real>& cons_forecast_state = (*forecast_state_at_lev)[IntVars::cons].array(mfi);
669  ApplyBndryForcing_Forecast(solverChoice, geom, tbx, tby, tbz, z_nd_arr,
670  xmom_src_arr, ymom_src_arr, zmom_src_arr,
671  rho_u, rho_v, rho_w,
672  rho_u_forecast_state, rho_v_forecast_state, rho_w_forecast_state,
673  cons_forecast_state);
674  }
675  if(solverChoice.init_type == InitType::HindCast and solverChoice.hindcast_surface_bcs) {
676  const Array4<const Real>& surface_state_arr = (*surface_state_at_lev).array(mfi);
678  xmom_src_arr, ymom_src_arr,
679  rho_u, rho_v,
680  cell_data, z_nd_arr,
681  surface_state_arr);
682  }
683  }
684 
685  // *****************************************************************************
686  // Real(8.) Add CANOPY source terms
687  // *****************************************************************************
688  if (solverChoice.do_forest_drag &&
689  ((is_slow_step && !use_canopy_fast) || (!is_slow_step && use_canopy_fast))) {
690  ParallelFor(tbx, tby, tbz,
691  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
692  {
693  const Real ux = u(i, j, k);
694  const Real uy = fourth * ( v(i, j , k ) + v(i-1, j , k )
695  + v(i, j+1, k ) + v(i-1, j+1, k ) );
696  const Real uz = fourth * ( w(i, j , k ) + w(i-1, j , k )
697  + w(i, j , k+1) + w(i-1, j , k+1) );
698  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
699  const Real f_drag = myhalf * (f_drag_arr(i, j, k) + f_drag_arr(i-1, j, k));
700  xmom_src_arr(i, j, k) -= f_drag * ux * windspeed;
701  },
702  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
703  {
704  const Real ux = fourth * ( u(i , j , k ) + u(i , j-1, k )
705  + u(i+1, j , k ) + u(i+1, j-1, k ) );
706  const Real uy = v(i, j, k);
707  const Real uz = fourth * ( w(i , j , k ) + w(i , j-1, k )
708  + w(i , j , k+1) + w(i , j-1, k+1) );
709  const amrex::Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
710  const Real f_drag = myhalf * (f_drag_arr(i, j, k) + f_drag_arr(i, j-1, k));
711  ymom_src_arr(i, j, k) -= f_drag * uy * windspeed;
712  },
713  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
714  {
715  const amrex::Real ux = fourth * ( u(i , j , k ) + u(i+1, j , k )
716  + u(i , j , k-1) + u(i+1, j , k-1) );
717  const amrex::Real uy = fourth * ( v(i , j , k ) + v(i , j+1, k )
718  + v(i , j , k-1) + v(i , j+1, k-1) );
719  const amrex::Real uz = w(i, j, k);
720  const amrex::Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
721  const Real f_drag = myhalf * (f_drag_arr(i, j, k) + f_drag_arr(i, j, k-1));
722  zmom_src_arr(i, j, k) -= f_drag * uz * windspeed;
723  });
724  }
725  // *****************************************************************************
726  // 9a. Add immersed source terms for terrain
727  // *****************************************************************************
728  if (solverChoice.terrain_type == TerrainType::ImmersedForcing &&
729  ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast))) {
730  // geometric properties
731  const Real* dx_arr = geom.CellSize();
732  const Real dx_x = dx_arr[0];
733  const Real dx_y = dx_arr[1];
734  const Real dx_z = dx_arr[2];
735 
736  const Real alpha_m = solverChoice.if_Cd_momentum;
737  const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, one/three);
739  const Real U_s = one; // unit velocity scale
740 
741  // MOST parameters
742  similarity_funs sfuns;
743  const Real ggg = CONST_GRAV;
744  const Real kappa = KAPPA;
745  const Real z0 = solverChoice.if_z0;
746  const Real tflux_in = solverChoice.if_surf_temp_flux;
747  const Real Olen_in = solverChoice.if_Olen_in;
748  const bool l_use_most = solverChoice.if_use_most;
749 
750  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
751  {
752  const Real ux = u(i, j, k);
753  const Real uy = fourth * ( v(i, j , k ) + v(i-1, j , k )
754  + v(i, j+1, k ) + v(i-1, j+1, k ) );
755  const Real uz = fourth * ( w(i, j , k ) + w(i-1, j , k )
756  + w(i, j , k+1) + w(i-1, j , k+1) );
757  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
758  const Real t_blank = myhalf * (t_blank_arr(i, j, k) + t_blank_arr(i-1, j, k));
759  const Real t_blank_above = myhalf * (t_blank_arr(i, j, k+1) + t_blank_arr(i-1, j, k+1));
760  const Real CdM = std::min(drag_coefficient / (windspeed + tiny), Real(1000.0));
761  const Real rho_xface = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i-1,j,k,Rho_comp) );
762 
763  if ((t_blank > 0 && (t_blank_above == zero)) && l_use_most) { // force to MOST value
764  // calculate tangential velocity one cell above
765  const Real ux2r = u(i, j, k+1) ;
766  const Real uy2r = fourth * ( v(i, j , k+1) + v(i-1, j , k+1)
767  + v(i, j+1, k+1) + v(i-1, j+1, k+1) ) ;
768  const Real h_windspeed2r = std::sqrt(ux2r * ux2r + uy2r * uy2r);
769 
770  // MOST
771  const Real theta_xface = (myhalf * (cell_data(i,j,k ,RhoTheta_comp) + cell_data(i-1,j,k, RhoTheta_comp))) / rho_xface;
772  const Real rho_xface_below = myhalf * ( cell_data(i,j,k-1,Rho_comp) + cell_data(i-1,j,k-1,Rho_comp) );
773  const Real theta_xface_below = (myhalf * (cell_data(i,j,k-1,RhoTheta_comp) + cell_data(i-1,j,k-1, RhoTheta_comp))) / rho_xface_below;
774  const Real theta_surf = theta_xface_below;
775 
776  Real psi_m = zero;
777  Real psi_h = zero;
778  Real ustar = h_windspeed2r * kappa / (std::log(Real(1.5) * dx_z / z0) - psi_m); // calculated from bottom of cell. Maintains flexibility for different Vf values
779  Real tflux = (tflux_in != Real(1e-8)) ? tflux_in : -(theta_xface - theta_surf) * ustar * kappa / (std::log(Real(1.5) * dx_z / z0) - psi_h);
780  Real Olen = (Olen_in != Real(1e-8)) ? Olen_in : -ustar * ustar * ustar * theta_xface / (kappa * ggg * tflux + tiny);
781  Real zeta = Real(1.5) * dx_z / Olen;
782 
783  // similarity functions
784  psi_m = sfuns.calc_psi_m(zeta);
785  psi_h = sfuns.calc_psi_h(zeta);
786  ustar = h_windspeed2r * kappa / (std::log(Real(1.5) * dx_z / z0) - psi_m);
787 
788  // prevent some unphysical math
789  if (!(ustar > zero && !std::isnan(ustar))) { ustar = zero; }
790  if (!(ustar < two && !std::isnan(ustar))) { ustar = two; }
791  if (psi_m > std::log(myhalf * dx_z / z0)) { psi_m = std::log(myhalf * dx_z / z0); }
792 
793  // determine target velocity
794  const Real uTarget = ustar / kappa * (std::log(myhalf * dx_z / z0) - psi_m);
795  Real uxTarget = uTarget * ux2r / (tiny + h_windspeed2r);
796  const Real bc_forcing_x = -(uxTarget - ux); // BC forcing pushes nonrelative velocity toward target velocity
797  xmom_src_arr(i, j, k) -= (1-t_blank) * rho_xface * CdM * U_s * bc_forcing_x; // if Vf low, force more strongly to MOST. If high, less forcing.
798  } else {
799  xmom_src_arr(i, j, k) -= t_blank * rho_xface * CdM * ux * windspeed;
800  }
801  });
802  ParallelFor(tby, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
803  {
804  const Real ux = fourth * ( u(i , j , k ) + u(i , j-1, k )
805  + u(i+1, j , k ) + u(i+1, j-1, k ) );
806  const Real uy = v(i, j, k);
807  const Real uz = fourth * ( w(i , j , k ) + w(i , j-1, k )
808  + w(i , j , k+1) + w(i , j-1, k+1) );
809  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
810  const Real t_blank = myhalf * (t_blank_arr(i, j, k) + t_blank_arr(i, j-1, k));
811  const Real t_blank_above = myhalf * (t_blank_arr(i, j, k+1) + t_blank_arr(i, j-1, k+1));
812  const Real CdM = std::min(drag_coefficient / (windspeed + tiny), Real(1000.0));
813  const Real rho_yface = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j-1,k,Rho_comp) );
814 
815  if ((t_blank > 0 && (t_blank_above == zero)) && l_use_most) { // force to MOST value
816  // calculate tangential velocity one cell above
817  const Real ux2r = fourth * ( u(i , j , k+1) + u(i , j-1, k+1)
818  + u(i+1, j , k+1) + u(i+1, j-1, k+1) );
819  const Real uy2r = v(i, j, k+1) ;
820  const Real h_windspeed2r = std::sqrt(ux2r * ux2r + uy2r * uy2r);
821 
822  // MOST
823  const Real theta_yface = (myhalf * (cell_data(i,j,k ,RhoTheta_comp) + cell_data(i,j-1,k, RhoTheta_comp))) / rho_yface;
824  const Real rho_yface_below = myhalf * ( cell_data(i,j,k-1,Rho_comp) + cell_data(i,j-1,k-1,Rho_comp) );
825  const Real theta_yface_below = (myhalf * (cell_data(i,j,k-1,RhoTheta_comp) + cell_data(i,j-1,k-1, RhoTheta_comp))) / rho_yface_below;
826  const Real theta_surf = theta_yface_below;
827 
828  Real psi_m = zero;
829  Real psi_h = zero;
830  Real ustar = h_windspeed2r * kappa / (std::log(Real(1.5) * dx_z / z0) - psi_m); // calculated from bottom of cell. Maintains flexibility for different Vf values
831  Real tflux = (tflux_in != Real(1e-8)) ? tflux_in : -(theta_yface - theta_surf) * ustar * kappa / (std::log(Real(1.5) * dx_z / z0) - psi_h);
832  Real Olen = (Olen_in != Real(1e-8)) ? Olen_in : -ustar * ustar * ustar * theta_yface / (kappa * ggg * tflux + tiny);
833  Real zeta = Real(1.5) * dx_z / Olen;
834 
835  // similarity functions
836  psi_m = sfuns.calc_psi_m(zeta);
837  psi_h = sfuns.calc_psi_h(zeta);
838  ustar = h_windspeed2r * kappa / (std::log(Real(1.5) * dx_z / z0) - psi_m);
839 
840  // prevent some unphysical math
841  if (!(ustar > zero && !std::isnan(ustar))) { ustar = zero; }
842  if (!(ustar < two && !std::isnan(ustar))) { ustar = two; }
843  if (psi_m > std::log(myhalf * dx_z / z0)) { psi_m = std::log(myhalf * dx_z / z0); }
844 
845  // determine target velocity
846  const Real uTarget = ustar / kappa * (std::log(myhalf * dx_z / z0) - psi_m);
847  Real uyTarget = uTarget * uy2r / (tiny + h_windspeed2r);
848  const Real bc_forcing_y = -(uyTarget - uy); // BC forcing pushes nonrelative velocity toward target velocity
849  ymom_src_arr(i, j, k) -= (1 - t_blank) * rho_yface * CdM * U_s * bc_forcing_y; // if Vf low, force more strongly to MOST. If high, less forcing.
850  } else {
851  ymom_src_arr(i, j, k) -= t_blank * rho_yface * CdM * uy * windspeed;
852  }
853  });
854  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
855  {
856  const Real ux = fourth * ( u(i , j , k ) + u(i+1, j , k )
857  + u(i , j , k-1) + u(i+1, j , k-1) );
858  const Real uy = fourth * ( v(i , j , k ) + v(i , j+1, k )
859  + v(i , j , k-1) + v(i , j+1, k-1) );
860  const Real uz = w(i, j, k);
861  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
862  const Real t_blank = myhalf * (t_blank_arr(i, j, k) + t_blank_arr(i, j, k-1));
863  const Real CdM = std::min(drag_coefficient / (windspeed + tiny), Real(1000.0));
864  const Real rho_zface = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp) );
865  zmom_src_arr(i, j, k) -= t_blank * rho_zface * CdM * uz * windspeed;
866  });
867  }
868 
869  // *****************************************************************************
870  // 9b. Add immersed source terms for buildings
871  // *****************************************************************************
872  if ((solverChoice.buildings_type == BuildingsType::ImmersedForcing ) &&
873  ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
874  {
875  // geometric properties
876  const Real* dx_arr = geom.CellSize();
877  const Real dx_x = dx_arr[0];
878  const Real dx_y = dx_arr[1];
879 
880  const Real alpha_m = solverChoice.if_Cd_momentum;
882  const Real min_t_blank = Real(0.005); // threshold for where immersed forcing acts
883 
884  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
885  {
886  const Real ux = u(i, j, k);
887  const Real uy = fourth * ( v(i, j , k ) + v(i-1, j , k )
888  + v(i, j+1, k ) + v(i-1, j+1, k ) );
889  const Real uz = fourth * ( w(i, j , k ) + w(i-1, j , k )
890  + w(i, j , k+1) + w(i-1, j , k+1) );
891  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
892 
893  Real t_blank = myhalf * (t_blank_arr(i, j, k) + t_blank_arr(i-1, j, k));
894  if (t_blank < min_t_blank) { t_blank = zero; }
895  const Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) : dx_arr[2];
896  const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, one/three);
897  const Real CdM = std::min(drag_coefficient / (windspeed + tiny), Real(1000.0));
898  const Real rho_xface = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i-1,j,k,Rho_comp) );
899  xmom_src_arr(i, j, k) -= t_blank * rho_xface * CdM * ux * windspeed;
900  });
901  ParallelFor(tby, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
902  {
903  const Real ux = fourth * ( u(i , j , k ) + u(i , j-1, k )
904  + u(i+1, j , k ) + u(i+1, j-1, k ) );
905  const Real uy = v(i, j, k);
906  const Real uz = fourth * ( w(i , j , k ) + w(i , j-1, k )
907  + w(i , j , k+1) + w(i , j-1, k+1) );
908  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
909 
910  Real t_blank = myhalf * (t_blank_arr(i, j, k) + t_blank_arr(i, j-1, k));
911  if (t_blank < min_t_blank) { t_blank = zero; }
912  const Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) : dx_arr[2];
913  const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, one/three);
914  const Real CdM = std::min(drag_coefficient / (windspeed + tiny), Real(1000.0));
915  const Real rho_yface = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j-1,k,Rho_comp) );
916  ymom_src_arr(i, j, k) -= t_blank * rho_yface * CdM * uy * windspeed;
917  });
918  ParallelFor(tbz, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
919  {
920  const Real ux = fourth * ( u(i , j , k ) + u(i+1, j , k )
921  + u(i , j , k-1) + u(i+1, j , k-1) );
922  const Real uy = fourth * ( v(i , j , k ) + v(i , j+1, k )
923  + v(i , j , k-1) + v(i , j+1, k-1) );
924  const Real uz = w(i, j, k);
925  const Real windspeed = std::sqrt(ux * ux + uy * uy + uz * uz);
926 
927  Real t_blank = myhalf * (t_blank_arr(i, j, k) + t_blank_arr(i, j, k-1));
928  if (t_blank < min_t_blank) { t_blank = zero; }
929  const Real dx_z = (z_nd_arr) ? (z_nd_arr(i,j,k) - z_nd_arr(i,j,k-1)) : dx_arr[2]; // ASW double check
930  const Real drag_coefficient = alpha_m / std::pow(dx_x*dx_y*dx_z, one/three);
931  const Real CdM = std::min(drag_coefficient / (windspeed + tiny), Real(1000.0));
932  const Real rho_zface = myhalf * ( cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp) );
933  zmom_src_arr(i, j, k) -= t_blank * rho_zface * CdM * uz * windspeed;
934  });
935  }
936 
937  // *****************************************************************************
938  // Real(10.) Enforce constant mass flux
939  // *****************************************************************************
940  if (is_slow_step && (enforce_massflux_x || enforce_massflux_y)) {
941  Real tau_inv = one / solverChoice.const_massflux_tau;
942 
943  ParallelFor(tbx, tby,
944  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
945  xmom_src_arr(i, j, k) += tau_inv * (rhoUA_target - rhoUA);
946  },
947  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
948  ymom_src_arr(i, j, k) += tau_inv * (rhoVA_target - rhoVA);
949  });
950  }
951 
952  } // mfi
953 }
void ApplyBndryForcing_Forecast(const SolverChoice &solverChoice, const Geometry geom, const Box &tbx, const Box &tby, const Box &tbz, const Array4< const Real > &z_phys_nd, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< Real > &rho_w_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &rho_w, const Array4< const Real > &rho_u_initial_state, const Array4< const Real > &rho_v_initial_state, const Array4< const Real > &rho_w_initial_state, const Array4< const Real > &cons_initial_state)
Definition: ERF_ApplyBndryForcing_Forecast.cpp:8
void ApplySpongeZoneBCsForMom(const SpongeChoice &spongeChoice, const Geometry geom, const Box &tbx, const Box &tby, const Box &tbz, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< Real > &rho_w_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &rho_w, const Array4< const Real > &r0, const Array4< const Real > &z_phys_nd, const Array4< const Real > &z_phys_cc)
Definition: ERF_ApplySpongeZoneBCs.cpp:169
void ApplySpongeZoneBCsForMom_ReadFromFile(const SpongeChoice &spongeChoice, const Geometry geom, const Box &tbx, const Box &tby, const Array4< const Real > &cell_data, const Array4< const Real > &z_phys_cc, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Vector< Real * > d_sponge_ptrs_at_lev)
Definition: ERF_ApplySpongeZoneBCs_ReadFromFile.cpp:8
void ApplySurfaceTreatment_BulkCoeff_Mom(const Box &tbx, const Box &tby, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &cons_state, const Array4< const Real > &z_phys_nd, const Array4< const Real > &surface_state_arr)
Definition: ERF_ApplySurfaceTreatment_BulkCoeff.cpp:8
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:39
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real PI
Definition: ERF_Constants.H:24
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:40
@ ubar
Definition: ERF_DataStruct.H:98
@ wbar
Definition: ERF_DataStruct.H:98
@ vbar
Definition: ERF_DataStruct.H:98
DirectionSelector< 2 > ZDir
Definition: ERF_DirectionSelector.H:38
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
Real z0
Definition: ERF_InitCustomPertVels_ScalarAdvDiff.H:8
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
void NumericalDiffusion_Ymom(const Box &bx, const Real dt, const Real num_diff_coeff, const Array4< const Real > &prim_data, const Array4< const Real > &cell_data, const Array4< Real > &rhs, const Array4< const Real > &mfx_arr, const Array4< const Real > &mfy_arr)
Definition: ERF_NumericalDiffusion.cpp:151
void NumericalDiffusion_Xmom(const Box &bx, const Real dt, const Real num_diff_coeff, const Array4< const Real > &prim_data, const Array4< const Real > &cell_data, const Array4< Real > &rhs, const Array4< const Real > &mfx_arr, const Array4< const Real > &mfy_arr)
Definition: ERF_NumericalDiffusion.cpp:85
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_PlaneAverage.H:14
eb_aux_ const * get_w_const_factory() const noexcept
Definition: ERF_EB.H:52
eb_aux_ const * get_v_const_factory() const noexcept
Definition: ERF_EB.H:51
eb_aux_ const * get_u_const_factory() const noexcept
Definition: ERF_EB.H:50
@ r0_comp
Definition: ERF_IndexDefines.H:73
@ ymom
Definition: ERF_IndexDefines.H:194
@ cons
Definition: ERF_IndexDefines.H:192
@ zmom
Definition: ERF_IndexDefines.H:195
@ xmom
Definition: ERF_IndexDefines.H:193
@ nr
Definition: ERF_Morrison.H:45
@ xvel
Definition: ERF_IndexDefines.H:175
@ cons
Definition: ERF_IndexDefines.H:174
@ yvel
Definition: ERF_IndexDefines.H:176
@ ww
Definition: ERF_AdvanceWSM6.cpp:105
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
real(kind=kind_phys), parameter, private r0
Definition: ERF_module_mp_wsm6.F90:21
bool rayleigh_damp_V
Definition: ERF_DampingStruct.H:85
amrex::Real rayleigh_dampcoef
Definition: ERF_DampingStruct.H:88
bool rayleigh_damp_W
Definition: ERF_DampingStruct.H:86
RayleighDampingType rayleigh_damping_type
Definition: ERF_DampingStruct.H:101
bool rayleigh_damp_U
Definition: ERF_DampingStruct.H:84
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > V_inp_sound_d
Definition: ERF_InputSoundingData.H:407
amrex::Vector< amrex::Real > input_sounding_time
Definition: ERF_InputSoundingData.H:394
amrex::Real tau_nudging
Definition: ERF_InputSoundingData.H:391
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > U_inp_sound_d
Definition: ERF_InputSoundingData.H:407
bool do_mom_advection
Definition: ERF_DataStruct.H:1239
amrex::Real coriolis_factor
Definition: ERF_DataStruct.H:1230
static MeshType mesh_type
Definition: ERF_DataStruct.H:1134
amrex::Real if_surf_temp_flux
Definition: ERF_DataStruct.H:1207
bool if_use_most
Definition: ERF_DataStruct.H:1211
DampingChoice dampingChoice
Definition: ERF_DataStruct.H:1144
amrex::Real const_massflux_v
Definition: ERF_DataStruct.H:1335
amrex::Real if_z0
Definition: ERF_DataStruct.H:1206
amrex::Real cosphi
Definition: ERF_DataStruct.H:1231
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > abl_geo_forcing
Definition: ERF_DataStruct.H:1307
bool hindcast_lateral_forcing
Definition: ERF_DataStruct.H:1344
int massflux_klo
Definition: ERF_DataStruct.H:1339
bool custom_w_subsidence
Definition: ERF_DataStruct.H:1237
bool nudging_from_input_sounding
Definition: ERF_DataStruct.H:1247
bool immersed_forcing_substep
Definition: ERF_DataStruct.H:1199
amrex::Real sinphi
Definition: ERF_DataStruct.H:1232
bool have_geo_wind_profile
Definition: ERF_DataStruct.H:1309
amrex::Real const_massflux_u
Definition: ERF_DataStruct.H:1334
amrex::Real if_Olen_in
Definition: ERF_DataStruct.H:1210
bool use_coriolis
Definition: ERF_DataStruct.H:1193
amrex::Real num_diff_coeff
Definition: ERF_DataStruct.H:1296
bool variable_coriolis
Definition: ERF_DataStruct.H:1311
amrex::Real if_Cd_momentum
Definition: ERF_DataStruct.H:1204
bool custom_forcing_prim_vars
Definition: ERF_DataStruct.H:1241
static BuildingsType buildings_type
Definition: ERF_DataStruct.H:1128
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1125
SpongeChoice spongeChoice
Definition: ERF_DataStruct.H:1145
static InitType init_type
Definition: ERF_DataStruct.H:1119
bool hindcast_surface_bcs
Definition: ERF_DataStruct.H:1345
bool do_forest_drag
Definition: ERF_DataStruct.H:1331
amrex::Real const_massflux_tau
Definition: ERF_DataStruct.H:1336
int massflux_khi
Definition: ERF_DataStruct.H:1340
bool forest_substep
Definition: ERF_DataStruct.H:1200
int ave_plane
Definition: ERF_DataStruct.H:1313
static SpongeType sponge_type
Definition: ERF_SpongeStruct.H:81
Definition: ERF_MOSTStress.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_m(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:90
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real calc_psi_h(amrex::Real zeta) const
Definition: ERF_MOSTStress.H:104
Here is the call graph for this function: