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