ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_MakeSources.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_SrcHeaders.H>
#include <ERF_TI_slow_headers.H>
#include <ERF_MOSTStress.H>
Include dependency graph for ERF_MakeSources.cpp:

Functions

void make_sources (int level, int, Real dt, Real time, const Vector< MultiFab > &S_data, const MultiFab &S_prim, MultiFab &source, const MultiFab &base_state, const MultiFab *z_phys_cc, const MultiFab &xvel, const MultiFab &yvel, const MultiFab *qheating_rates, MultiFab *terrain_blank, const Geometry geom, const SolverChoice &solverChoice, Vector< std::unique_ptr< MultiFab >> &mapfac, const MultiFab *rhotheta_src, const MultiFab *rhoqt_src, const Real *dptr_wbar_sub, const Vector< Real * > d_rayleigh_ptrs_at_lev, const Real *d_sinesq_at_lev, const MultiFab *surface_state_at_lev, InputSoundingData &input_sounding_data, TurbulentPerturbation &turbPert, bool is_slow_step)
 

Function Documentation

◆ make_sources()

void make_sources ( int  level,
int  ,
Real  dt,
Real  time,
const Vector< MultiFab > &  S_data,
const MultiFab &  S_prim,
MultiFab &  source,
const MultiFab &  base_state,
const MultiFab *  z_phys_cc,
const MultiFab &  xvel,
const MultiFab &  yvel,
const MultiFab *  qheating_rates,
MultiFab *  terrain_blank,
const Geometry  geom,
const SolverChoice solverChoice,
Vector< std::unique_ptr< MultiFab >> &  mapfac,
const MultiFab *  rhotheta_src,
const MultiFab *  rhoqt_src,
const Real dptr_wbar_sub,
const Vector< Real * >  d_rayleigh_ptrs_at_lev,
const Real d_sinesq_at_lev,
const MultiFab *  surface_state_at_lev,
InputSoundingData input_sounding_data,
TurbulentPerturbation turbPert,
bool  is_slow_step 
)

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

Parameters
[in]levellevel of resolution
[in]nrkwhich RK stage
[in]dtslow time step
[in]S_datacurrent solution
[in]S_primprimitive variables (i.e. conserved variables divided by density)
[in]sourcesource terms for conserved variables
[in]geomContainer for geometric information
[in]solverChoiceContainer for solver parameters
[in]mapfacmap factors
[in]dptr_rhotheta_srccustom temperature source term
[in]dptr_rhoqt_srccustom moisture source term
[in]dptr_wbar_subsubsidence source term
[in]d_rayleigh_ptrs_at_levVector of {strength of Rayleigh damping, reference value of theta} used to define Rayleigh damping
[in]d_sinesq_at_levsin( (pi/2) (z-z_t)/(damping depth)) at cell centers
58 {
59  BL_PROFILE_REGION("erf_make_sources()");
60 
61  // *****************************************************************************
62  // Initialize source to zero since we re-compute it every RK stage
63  // *****************************************************************************
64  if (is_slow_step) {
65  source.setVal(0.);
66  } else {
67  source.setVal(0.0,Rho_comp,2);
68  }
69 
70  const bool l_use_ndiff = solverChoice.use_num_diff;
71 
72  TurbChoice tc = solverChoice.turbChoice[level];
73  const bool l_use_KE = tc.use_tke;
74  const bool l_diff_KE = tc.diffuse_tke_3D;
75 
76  const Box& domain = geom.Domain();
77 
78  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
79  const GpuArray<Real, AMREX_SPACEDIM> dx = geom.CellSizeArray();
80 
81  MultiFab r_hse (base_state, make_alias, BaseState::r0_comp , 1);
82 
83  Real* thetabar = d_rayleigh_ptrs_at_lev[Rayleigh::thetabar];
84 
85  // flags to apply certain source terms in substep call only
86  bool use_Rayleigh_fast = ( (solverChoice.dampingChoice.rayleigh_damping_type == RayleighDampingType::FastExplicit) ||
88  bool use_ImmersedForcing_fast = solverChoice.immersed_forcing_substep;
89 
90  // flag for a moisture model
91  bool has_moisture = (solverChoice.moisture_type != MoistureType::None);
92 
93  // *****************************************************************************
94  // Planar averages for subsidence terms
95  // *****************************************************************************
96  Table1D<Real> dptr_r_plane, dptr_t_plane, dptr_qv_plane, dptr_qc_plane;
97  TableData<Real, 1> r_plane_tab, t_plane_tab, qv_plane_tab, qc_plane_tab;
98  bool compute_averages = false;
99  compute_averages = compute_averages ||
100  ( is_slow_step && (dptr_wbar_sub || solverChoice.nudging_from_input_sounding) );
101 
102  if (compute_averages)
103  {
104  //
105  // The call to "compute_averages" currently does all the components in one call
106  // We can then extract each component separately with the "line_average" call
107  //
108  // We need just one ghost cell in the vertical
109  //
110  IntVect ng_c(S_data[IntVars::cons].nGrowVect()); ng_c[2] = 1;
111  //
112  // With no moisture we only (rho) and (rho theta); with moisture we also do qv and qc
113  // We use the alias here to control ncomp inside the PlaneAverage
114  //
115  int ncomp = (!has_moisture) ? 2 : RhoQ2_comp+1;
116  MultiFab cons(S_data[IntVars::cons], make_alias, 0, ncomp);
117 
118  PlaneAverage cons_ave(&cons, geom, solverChoice.ave_plane, ng_c);
119  cons_ave.compute_averages(ZDir(), cons_ave.field());
120 
121  int ncell = cons_ave.ncell_line();
122 
123  Gpu::HostVector< Real> r_plane_h(ncell);
124  Gpu::DeviceVector< Real> r_plane_d(ncell);
125 
126  Gpu::HostVector< Real> t_plane_h(ncell);
127  Gpu::DeviceVector< Real> t_plane_d(ncell);
128 
129  cons_ave.line_average(Rho_comp , r_plane_h);
130  cons_ave.line_average(RhoTheta_comp, t_plane_h);
131 
132  Gpu::copyAsync(Gpu::hostToDevice, r_plane_h.begin(), r_plane_h.end(), r_plane_d.begin());
133  Gpu::copyAsync(Gpu::hostToDevice, t_plane_h.begin(), t_plane_h.end(), t_plane_d.begin());
134 
135  Real* dptr_r = r_plane_d.data();
136  Real* dptr_t = t_plane_d.data();
137 
138  Box tdomain = domain; tdomain.grow(2,ng_c[2]);
139  r_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
140  t_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
141 
142  int offset = ng_c[2];
143 
144  dptr_r_plane = r_plane_tab.table();
145  dptr_t_plane = t_plane_tab.table();
146  ParallelFor(ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
147  {
148  dptr_r_plane(k-offset) = dptr_r[k];
149  dptr_t_plane(k-offset) = dptr_t[k];
150  });
151 
152  if (has_moisture)
153  {
154  Gpu::HostVector< Real> qv_plane_h(ncell), qc_plane_h(ncell);
155  Gpu::DeviceVector<Real> qv_plane_d(ncell), qc_plane_d(ncell);
156 
157  // Water vapor
158  cons_ave.line_average(RhoQ1_comp, qv_plane_h);
159  Gpu::copyAsync(Gpu::hostToDevice, qv_plane_h.begin(), qv_plane_h.end(), qv_plane_d.begin());
160 
161  // Cloud water
162  cons_ave.line_average(RhoQ2_comp, qc_plane_h);
163  Gpu::copyAsync(Gpu::hostToDevice, qc_plane_h.begin(), qc_plane_h.end(), qc_plane_d.begin());
164 
165  Real* dptr_qv = qv_plane_d.data();
166  Real* dptr_qc = qc_plane_d.data();
167 
168  qv_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
169  qc_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
170 
171  dptr_qv_plane = qv_plane_tab.table();
172  dptr_qc_plane = qc_plane_tab.table();
173  ParallelFor(ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
174  {
175  dptr_qv_plane(k-offset) = dptr_qv[k];
176  dptr_qc_plane(k-offset) = dptr_qc[k];
177  });
178  }
179  }
180 
181  // *****************************************************************************
182  // Radiation flux vector for four stream approximation
183  // *****************************************************************************
184  // NOTE: The fluxes live on w-faces
185  int klo = domain.smallEnd(0);
186  int khi = domain.bigEnd(2);
187  int nk = khi - klo + 2;
188  Gpu::DeviceVector<Real> radiation_flux(nk,0.0);
189  Gpu::DeviceVector<Real> q_integral(nk,0.0);
190  Real* rad_flux = radiation_flux.data();
191  Real* q_int = q_integral.data();
192 
193  // *****************************************************************************
194  // Define source term for cell-centered conserved variables, from
195  // 1. user-defined source terms for (rho theta) and (rho q_t)
196  // 2. radiation for (rho theta)
197  // 3. Rayleigh damping for (rho theta)
198  // 4. custom forcing for (rho theta) and (rho Q1)
199  // 5. custom subsidence for (rho theta) and (rho Q1)
200  // 6. numerical diffusion for (rho theta)
201  // 7. sponging
202  // 8. turbulent perturbation
203  // 9. nudging towards input sounding values (only for theta)
204  // 10a. Immersed forcing for terrain
205  // 10b. Immersed forcing for buildings
206  // 11. Four stream radiation source for (rho theta)
207  // *****************************************************************************
208 
209  // ***********************************************************************************************
210  // Add remaining source terms
211  // ***********************************************************************************************
212 #ifdef _OPENMP
213 #pragma omp parallel if (Gpu::notInLaunchRegion())
214 #endif
215  {
216  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
217  {
218  Box bx = mfi.tilebox();
219 
220  const Array4<const Real>& cell_data = S_data[IntVars::cons].array(mfi);
221  const Array4<const Real>& cell_prim = S_prim.array(mfi);
222  const Array4<Real> & cell_src = source.array(mfi);
223 
224  const Array4<const Real>& r0 = r_hse.const_array(mfi);
225 
226  const Array4<const Real>& z_cc_arr = z_phys_cc->const_array(mfi);
227 
228  const Array4<const Real>& t_blank_arr = (terrain_blank) ? terrain_blank->const_array(mfi) :
229  Array4<const Real>{};
230 
231 
232  // *************************************************************************************
233  // 2. Add radiation source terms to (rho theta)
234  // *************************************************************************************
235  if (solverChoice.rad_type != RadiationType::None && is_slow_step) {
236  auto const& qheating_arr = qheating_rates->const_array(mfi);
237  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
238  {
239  // Short-wavelength and long-wavelength radiation source terms
240  cell_src(i,j,k,RhoTheta_comp) += cell_data(i,j,k,Rho_comp) * ( qheating_arr(i,j,k,0) + qheating_arr(i,j,k,1) );
241  });
242  }
243 
244 
245  // *************************************************************************************
246  // 3. Add Rayleigh damping for (rho theta)
247  // *************************************************************************************
248  Real dampcoef = solverChoice.dampingChoice.rayleigh_dampcoef;
249 
250  if (solverChoice.dampingChoice.rayleigh_damp_T) {
251  if ((is_slow_step && !use_Rayleigh_fast) || (!is_slow_step && use_Rayleigh_fast)) {
252  int n = RhoTheta_comp;
253  int nr = Rho_comp;
254  int np = PrimTheta_comp;
255 
256  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
257  {
258  Real theta = cell_prim(i,j,k,np);
259  Real sinesq = d_sinesq_at_lev[k];
260  cell_src(i, j, k, n) -= dampcoef*sinesq * (theta - thetabar[k]) * cell_data(i,j,k,nr);
261  });
262  }
263  }
264 
265  // *************************************************************************************
266  // 4. Add custom forcing for (rho theta)
267  // *************************************************************************************
268  if (solverChoice.custom_rhotheta_forcing && is_slow_step) {
269  const int n = RhoTheta_comp;
270  auto const& rhotheta_src_arr = rhotheta_src->const_array(mfi);
271  if (solverChoice.spatial_rhotheta_forcing)
272  {
273  if (solverChoice.custom_forcing_prim_vars) {
274  const int nr = Rho_comp;
275  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
276  {
277  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhotheta_src_arr(i, j, k);
278  });
279  } else {
280  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
281  {
282  cell_src(i, j, k, n) += rhotheta_src_arr(i, j, k);
283  });
284  }
285  } else {
286  if (solverChoice.custom_forcing_prim_vars) {
287  const int nr = Rho_comp;
288  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
289  {
290  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhotheta_src_arr(0, 0, k);
291  });
292  } else {
293  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
294  {
295  cell_src(i, j, k, n) += rhotheta_src_arr(0, 0, k);
296  });
297  }
298  }
299  }
300 
301  // *************************************************************************************
302  // 4. Add custom forcing for RhoQ1
303  // *************************************************************************************
304  if (solverChoice.custom_moisture_forcing && is_slow_step) {
305  const int n = RhoQ1_comp;
306  auto const& rhoqt_src_arr = rhoqt_src->const_array(mfi);
307  if (solverChoice.spatial_moisture_forcing)
308  {
309  if (solverChoice.custom_forcing_prim_vars) {
310  const int nr = Rho_comp;
311  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
312  {
313  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhoqt_src_arr(i, j, k);
314  });
315  } else {
316  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
317  {
318  cell_src(i, j, k, n) += rhoqt_src_arr(i, j, k);
319  });
320  }
321  } else {
322  if (solverChoice.custom_forcing_prim_vars) {
323  const int nr = Rho_comp;
324  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
325  {
326  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhoqt_src_arr(0, 0, k);
327  });
328  } else {
329  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
330  {
331  cell_src(i, j, k, n) += rhoqt_src_arr(0, 0, k);
332  });
333  }
334  }
335  }
336 
337  // *************************************************************************************
338  // 5. Add custom subsidence for (rho theta)
339  // *************************************************************************************
340  if (solverChoice.custom_w_subsidence && is_slow_step && solverChoice.do_theta_advection) {
341  const int n = RhoTheta_comp;
342  if (solverChoice.custom_forcing_prim_vars) {
343  const int nr = Rho_comp;
344  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
345  {
346  Real dzInv = (z_cc_arr) ? 1.0/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : 0.5*dxInv[2];
347  Real T_hi = dptr_t_plane(k+1) / dptr_r_plane(k+1);
348  Real T_lo = dptr_t_plane(k-1) / dptr_r_plane(k-1);
349  Real wbar_cc = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
350  cell_src(i, j, k, n) -= cell_data(i,j,k,nr) * wbar_cc * (T_hi - T_lo) * dzInv;
351  });
352  } else {
353  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
354  {
355  Real dzInv = (z_cc_arr) ? 1.0/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : 0.5*dxInv[2];
356  Real T_hi = dptr_t_plane(k+1) / dptr_r_plane(k+1);
357  Real T_lo = dptr_t_plane(k-1) / dptr_r_plane(k-1);
358  Real wbar_cc = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
359  cell_src(i, j, k, n) -= wbar_cc * (T_hi - T_lo) * dzInv;
360  });
361  }
362  }
363 
364  // *************************************************************************************
365  // 5. Add custom subsidence for RhoQ1 and RhoQ2
366  // *************************************************************************************
367  if (solverChoice.custom_w_subsidence && (solverChoice.moisture_type != MoistureType::None) && is_slow_step) {
368  const int nv = RhoQ1_comp;
369  if (solverChoice.custom_forcing_prim_vars) {
370  const int nr = Rho_comp;
371  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
372  {
373  Real dzInv = (z_cc_arr) ? 1.0/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : 0.5*dxInv[2];
374  Real Qv_hi = dptr_qv_plane(k+1) / dptr_r_plane(k+1);
375  Real Qv_lo = dptr_qv_plane(k-1) / dptr_r_plane(k-1);
376  Real Qc_hi = dptr_qc_plane(k+1) / dptr_r_plane(k+1);
377  Real Qc_lo = dptr_qc_plane(k-1) / dptr_r_plane(k-1);
378  Real wbar_cc = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
379  cell_src(i, j, k, nv ) -= cell_data(i,j,k,nr) * wbar_cc * (Qv_hi - Qv_lo) * dzInv;
380  cell_src(i, j, k, nv+1) -= cell_data(i,j,k,nr) * wbar_cc * (Qc_hi - Qc_lo) * dzInv;
381  });
382  } else {
383  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
384  {
385  Real dzInv = (z_cc_arr) ? 1.0/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : 0.5*dxInv[2];
386  Real Qv_hi = dptr_qv_plane(k+1) / dptr_r_plane(k+1);
387  Real Qv_lo = dptr_qv_plane(k-1) / dptr_r_plane(k-1);
388  Real Qc_hi = dptr_qc_plane(k+1) / dptr_r_plane(k+1);
389  Real Qc_lo = dptr_qc_plane(k-1) / dptr_r_plane(k-1);
390  Real wbar_cc = 0.5 * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
391  cell_src(i, j, k, nv ) -= wbar_cc * (Qv_hi - Qv_lo) * dzInv;
392  cell_src(i, j, k, nv+1) -= wbar_cc * (Qc_hi - Qc_lo) * dzInv;
393  });
394  }
395  }
396 
397  // *************************************************************************************
398  // 6. Add numerical diffuion for rho and (rho theta)
399  // *************************************************************************************
400  if (l_use_ndiff && is_slow_step) {
401  int sc;
402  int nc;
403 
404  const Array4<const Real>& mf_mx = mapfac[MapFacType::m_x]->const_array(mfi);
405  const Array4<const Real>& mf_my = mapfac[MapFacType::m_y]->const_array(mfi);
406 
407  // Rho is a special case
408  NumericalDiffusion_Scal(bx, sc=0, nc=1, dt, solverChoice.num_diff_coeff,
409  cell_data, cell_data, cell_src, mf_mx, mf_my);
410 
411  // Other scalars proceed as normal
412  NumericalDiffusion_Scal(bx, sc=1, nc=1, dt, solverChoice.num_diff_coeff,
413  cell_prim, cell_data, cell_src, mf_mx, mf_my);
414 
415 
416  if (l_use_KE && l_diff_KE) {
417  NumericalDiffusion_Scal(bx, sc=RhoKE_comp, nc=1, dt, solverChoice.num_diff_coeff,
418  cell_prim, cell_data, cell_src, mf_mx, mf_my);
419  }
420 
422  cell_prim, cell_data, cell_src, mf_mx, mf_my);
423  }
424 
425  // *************************************************************************************
426  // 7. Add sponging
427  // *************************************************************************************
428  if(!(solverChoice.spongeChoice.sponge_type == "input_sponge") && is_slow_step){
429  ApplySpongeZoneBCsForCC(solverChoice.spongeChoice, geom, bx, cell_src, cell_data, r0, z_cc_arr);
430  }
431 
432  if(solverChoice.init_type == InitType::HindCast and solverChoice.hindcast_surface_bcs) {
433  const Array4<const Real>& surface_state_arr = (*surface_state_at_lev).array(mfi);
434  ApplySurfaceTreatment_BulkCoeff_CC(bx, cell_src, cell_data, z_cc_arr, surface_state_arr);
435  }
436 
437 
438  // *************************************************************************************
439  // 8. Add perturbation
440  // *************************************************************************************
441  if (solverChoice.pert_type == PerturbationType::Source && is_slow_step) {
442  auto m_ixtype = S_data[IntVars::cons].boxArray().ixType(); // Conserved term
443  const amrex::Array4<const amrex::Real>& pert_cell = turbPert.pb_cell[level].const_array(mfi);
444  turbPert.apply_tpi(level, bx, RhoTheta_comp, m_ixtype, cell_src, pert_cell); // Applied as source term
445  }
446 
447  // *************************************************************************************
448  // 9. Add nudging towards value specified in input sounding
449  // *************************************************************************************
450  if (solverChoice.nudging_from_input_sounding && is_slow_step)
451  {
452  int itime_n = 0;
453  int itime_np1 = 0;
454  Real coeff_n = Real(1.0);
455  Real coeff_np1 = Real(0.0);
456 
457  Real tau_inv = Real(1.0) / input_sounding_data.tau_nudging;
458 
459  int n_sounding_times = input_sounding_data.input_sounding_time.size();
460 
461  for (int nt = 1; nt < n_sounding_times; nt++) {
462  if (time > input_sounding_data.input_sounding_time[nt]) itime_n = nt;
463  }
464  if (itime_n == n_sounding_times-1) {
465  itime_np1 = itime_n;
466  } else {
467  itime_np1 = itime_n+1;
468  coeff_np1 = (time - input_sounding_data.input_sounding_time[itime_n]) /
469  (input_sounding_data.input_sounding_time[itime_np1] - input_sounding_data.input_sounding_time[itime_n]);
470  coeff_n = Real(1.0) - coeff_np1;
471  }
472 
473  const Real* theta_inp_sound_n = input_sounding_data.theta_inp_sound_d[itime_n].dataPtr();
474  const Real* theta_inp_sound_np1 = input_sounding_data.theta_inp_sound_d[itime_np1].dataPtr();
475 
476  const int n = RhoTheta_comp;
477  const int nr = Rho_comp;
478 
479  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
480  {
481  Real nudge = (coeff_n*theta_inp_sound_n[k] + coeff_np1*theta_inp_sound_np1[k]) - (dptr_t_plane(k)/dptr_r_plane(k));
482  nudge *= tau_inv;
483  cell_src(i, j, k, n) += cell_data(i, j, k, nr) * nudge;
484  });
485  }
486 
487  // *************************************************************************************
488  // 10a. Add immersed source terms for terrain
489  // *************************************************************************************
490  if (solverChoice.terrain_type == TerrainType::ImmersedForcing &&
491  ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
492  {
493  const Array4<const Real>& u = xvel.array(mfi);
494  const Array4<const Real>& v = yvel.array(mfi);
495 
496  // geometric properties
497  const Real* dx_arr = geom.CellSize();
498  const Real dx_x = dx_arr[0];
499  const Real dx_y = dx_arr[1];
500  const Real dx_z = dx_arr[2];
501 
502  const Real alpha_h = solverChoice.if_Cd_scalar;
503  const Real drag_coefficient = alpha_h / std::pow(dx_x*dx_y*dx_z, 1./3.);
505  const Real U_s = 1.0; // unit velocity scale
506 
507  // MOST parameters
508  similarity_funs sfuns;
509  const Real ggg = CONST_GRAV;
510  const Real kappa = KAPPA;
511  const Real z0 = solverChoice.if_z0;
512  const Real tflux = solverChoice.if_surf_temp_flux;
513  const Real init_surf_temp = solverChoice.if_init_surf_temp;
514  const Real surf_heating_rate = solverChoice.if_surf_heating_rate;
515  const Real Olen_in = solverChoice.if_Olen_in;
516 
517  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
518  {
519  const Real t_blank = t_blank_arr(i, j, k);
520  const Real t_blank_above = t_blank_arr(i, j, k+1);
521  const Real ux_cc_2r = 0.5 * (u(i ,j ,k+1) + u(i+1,j ,k+1));
522  const Real uy_cc_2r = 0.5 * (v(i ,j ,k+1) + v(i ,j+1,k+1));
523  const Real h_windspeed2r = std::sqrt(ux_cc_2r * ux_cc_2r + uy_cc_2r * uy_cc_2r);
524 
525  const Real theta = cell_data(i,j,k ,RhoTheta_comp) / cell_data(i,j,k ,Rho_comp);
526  const Real theta_neighbor = cell_data(i,j,k+1,RhoTheta_comp) / cell_data(i,j,k+1,Rho_comp);
527 
528  // SURFACE TEMP AND HEATING/COOLING RATE
529  if (init_surf_temp > 0.0) {
530  if (t_blank > 0 && (t_blank_above == 0.0)) { // force to MOST value
531  const Real surf_temp = init_surf_temp + surf_heating_rate*time/3600;
532  const Real bc_forcing_rt_srf = -(cell_data(i,j,k-1,Rho_comp) * surf_temp - cell_data(i,j,k-1,RhoTheta_comp));
533  cell_src(i, j, k-1, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf; // k-1
534  }
535  }
536 
537  // SURFACE HEAT FLUX
538  if (tflux != 1e-8){
539  if (t_blank > 0 && (t_blank_above == 0.0)) { // force to MOST value
540  Real psi_m = 0.0;
541  Real psi_h = 0.0;
542  Real psi_h_neighbor = 0.0;
543  Real ustar = h_windspeed2r * kappa / (std::log((1.5) * dx_z / z0) - psi_m);
544  const Real Olen = -ustar * ustar * ustar * theta / (kappa * ggg * tflux + tiny);
545  const Real zeta = (0.5) * dx_z / Olen;
546  const Real zeta_neighbor = (1.5) * dx_z / Olen;
547 
548  // similarity functions
549  psi_m = sfuns.calc_psi_m(zeta);
550  psi_h = sfuns.calc_psi_h(zeta);
551  psi_h_neighbor = sfuns.calc_psi_h(zeta_neighbor);
552  ustar = h_windspeed2r * kappa / (std::log((1.5) * dx_z / z0) - psi_m);
553 
554  // prevent some unphysical math
555  if (!(ustar > 0.0 && !std::isnan(ustar))) { ustar = 0.0; }
556  if (!(ustar < 2.0 && !std::isnan(ustar))) { ustar = 2.0; }
557  if (psi_h_neighbor > std::log(1.5 * dx_z / z0)) { psi_h_neighbor = std::log(1.5 * dx_z / z0); }
558  if (psi_h > std::log(0.5 * dx_z / z0)) { psi_h = std::log(0.5 * dx_z / z0); }
559 
560  // We do not know the actual temperature so use cell above
561  const Real thetastar = theta * ustar * ustar / (kappa * ggg * Olen);
562  const Real surf_temp = theta_neighbor - thetastar / kappa * (std::log((1.5) * dx_z / z0) - psi_h_neighbor);
563  const Real tTarget = surf_temp + thetastar / kappa * (std::log((0.5) * dx_z / z0) - psi_h);
564 
565  const Real bc_forcing_rt = -(cell_data(i,j,k,Rho_comp) * tTarget - cell_data(i,j,k,RhoTheta_comp));
566  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt;
567  }
568  }
569 
570  // OBUKHOV LENGTH
571  if (Olen_in != 1e-8){
572  if (t_blank > 0 && (t_blank_above == 0.0)) { // force to MOST value
573  const Real Olen = Olen_in;
574  const Real zeta = (0.5) * dx_z / Olen;
575  const Real zeta_neighbor = (1.5) * dx_z / Olen;
576 
577  // similarity functions
578  const Real psi_m = sfuns.calc_psi_m(zeta);
579  const Real psi_h = sfuns.calc_psi_h(zeta);
580  const Real psi_h_neighbor = sfuns.calc_psi_h(zeta_neighbor);
581  const Real ustar = h_windspeed2r * kappa / (std::log((1.5) * dx_z / z0) - psi_m);
582 
583  // We do not know the actual temperature so use cell above
584  const Real thetastar = theta * ustar * ustar / (kappa * ggg * Olen);
585  const Real surf_temp = theta_neighbor - thetastar / kappa * (std::log((1.5) * dx_z / z0) - psi_h_neighbor);
586  const Real tTarget = surf_temp + thetastar / kappa * (std::log((0.5) * dx_z / z0) - psi_h);
587 
588  const Real bc_forcing_rt = -(cell_data(i,j,k,Rho_comp) * tTarget - cell_data(i,j,k,RhoTheta_comp));
589  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt;
590  }
591  }
592 
593  });
594  }
595 
596  // *************************************************************************************
597  // 10b. Add immersed source terms for buildings
598  // *************************************************************************************
599  if ((solverChoice.buildings_type == BuildingsType::ImmersedForcing ) &&
600  ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
601  {
602  // geometric properties
603  const Real* dx_arr = geom.CellSize();
604  const Real dx_x = dx_arr[0];
605  const Real dx_y = dx_arr[1];
606 
607  const Real alpha_h = solverChoice.if_Cd_scalar;
608  const Real U_s = 1.0; // unit velocity scale
609  const Real min_t_blank = 0.005;
610 
611  const Real init_surf_temp = solverChoice.if_init_surf_temp;
612  const Real surf_heating_rate = solverChoice.if_surf_heating_rate;
613 
614  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
615  {
616  Real t_blank = t_blank_arr(i, j, k);
617  Real t_blank_below = t_blank_arr(i, j, k-1);
618  Real t_blank_above = t_blank_arr(i, j, k+1);
619  Real t_blank_north = t_blank_arr(i , j+1, k);
620  Real t_blank_south = t_blank_arr(i , j-1, k);
621  Real t_blank_east = t_blank_arr(i+1, j , k);
622  Real t_blank_west = t_blank_arr(i-1, j , k);
623  if (t_blank < min_t_blank) { t_blank = 0.0; } // deal with situations where very small volfrac exist
624  if (t_blank_below < min_t_blank) { t_blank_below = 0.0; }
625  if (t_blank_north < min_t_blank) { t_blank_north = 0.0; }
626  if (t_blank_south < min_t_blank) { t_blank_south = 0.0; }
627  if (t_blank_east < min_t_blank) { t_blank_east = 0.0; }
628  if (t_blank_west < min_t_blank) { t_blank_west = 0.0; }
629 
630  Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) : dx[2];
631  Real drag_coefficient = alpha_h / std::pow(dx_x*dx_y*dx_z, 1./3.);
632 
633  // SURFACE TEMP AND HEATING/COOLING RATE
634  if (init_surf_temp > 0.0) {
635  const Real surf_temp = init_surf_temp + surf_heating_rate*time/3600;
636  if (t_blank > 0 && (t_blank_above == 0.0) && (t_blank_below == 1.0)) { // building roof
637  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
638  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
639 
640  } else if (((t_blank > 0 && t_blank < t_blank_west && t_blank_east == 0.0) ||
641  (t_blank > 0 && t_blank < t_blank_east && t_blank_west == 0.0) ||
642  (t_blank > 0 && t_blank < t_blank_north && t_blank_south == 0.0) ||
643  (t_blank > 0 && t_blank < t_blank_south && t_blank_north == 0.0))) {
644  // this should enter for just building walls
645  // walls are currently separated to allow for flexible in the future to heat walls differently
646 
647  // south face
648  if ((t_blank < t_blank_north) && (t_blank_north == 1.0)) {
649  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
650  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
651  }
652 
653  // north face
654  if ((t_blank < t_blank_south) && (t_blank_south == 1.0)) {
655  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
656  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
657  }
658 
659  // west face
660  if ((t_blank < t_blank_east) && (t_blank_east == 1.0)) {
661  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
662  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
663  }
664 
665  // east face
666  if ((t_blank < t_blank_west) && (t_blank_west == 1.0)) {
667  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
668  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
669  }
670 
671  }
672  }
673  });
674  }
675 
676  // *************************************************************************************
677  // 11. Add 4 stream radiation src to RhoTheta
678  // *************************************************************************************
679  if (solverChoice.four_stream_radiation && has_moisture && is_slow_step)
680  {
681  AMREX_ALWAYS_ASSERT((bx.smallEnd(2) == klo) && (bx.bigEnd(2) == khi));
682  Real D = 3.75e-6; // [s^-1]
683  Real F0 = 70; // [W/m^2]
684  Real F1 = 22; // [W/m^2]
685  Real krad = 85; // [m^2 kg^-1]
686  Real qt_i = 0.008;
687 
688  Box xybx = makeSlab(bx,2,klo);
689  ParallelFor(xybx, [=] AMREX_GPU_DEVICE(int i, int j, int /*k*/) noexcept
690  {
691  // Inclusive scan at w-faces for the Q integral (also find "i" values)
692  q_int[0] = 0.0;
693  Real zi = 0.5 * (z_cc_arr(i,j,khi) + z_cc_arr(i,j,khi-1));
694  Real rhoi = 0.5 * (cell_data(i,j,khi,Rho_comp) + cell_data(i,j,khi-1,Rho_comp));
695  for (int k(klo+1); k<=khi+1; ++k) {
696  int lk = k - klo;
697  // Average to w-faces when looping w-faces
698  Real dz = (z_cc_arr) ? 0.5 * (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-2)) : dx[2];
699  q_int[lk] = q_int[lk-1] + krad * cell_data(i,j,k-1,Rho_comp) * cell_data(i,j,k-1,RhoQ2_comp) * dz;
700  Real qt_hi = cell_data(i,j,k ,RhoQ1_comp) + cell_data(i,j,k ,RhoQ2_comp);
701  Real qt_lo = cell_data(i,j,k-1,RhoQ1_comp) + cell_data(i,j,k-1,RhoQ2_comp);
702  if ( (qt_lo > qt_i) && (qt_hi < qt_i) ) {
703  zi = 0.5 * (z_cc_arr(i,j,k) + z_cc_arr(i,j,k-1));
704  rhoi = 0.5 * (cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp));
705  }
706  }
707 
708  // Decompose the integral to get the fluxes at w-faces
709  Real q_int_inf = q_int[khi+1];
710  for (int k(klo); k<=khi+1; ++k) {
711  int lk = k - klo;
712  Real z = 0.5 * (z_cc_arr(i,j,k) + z_cc_arr(i,j,k-1));
713  rad_flux[lk] = F1*std::exp(-q_int[lk]) + F0*std::exp(-(q_int_inf - q_int[lk]));
714  if (z > zi) {
715  rad_flux[lk] += rhoi * Cp_d * D * ( std::pow(z-zi,4./3.)/4. + zi*std::pow(z-zi,1./3.) ) ;
716  }
717  }
718 
719  // Compute the radiative heating source
720  for (int k(klo); k<=khi; ++k) {
721  int lk = k - klo;
722  // Average to w-faces when looping CC
723  Real dzInv = (z_cc_arr) ? 1.0/ (0.5 * (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1))) : dxInv[2];
724  // NOTE: Fnet = Up - Dn (all fluxes are up here)
725  // dT/dt = dF/dz * (1/(-rho*Cp))
726  Real dTdt = (rad_flux[lk+1] - rad_flux[lk]) * dzInv / (-cell_data(i,j,k,Rho_comp)*Cp_d);
727  Real qv = cell_data(i,j,k,RhoQ1_comp)/cell_data(i,j,k,Rho_comp);
728  Real iexner = 1./getExnergivenRTh(cell_data(i,j,k,RhoTheta_comp), R_d/Cp_d, qv);
729  // Convert dT/dt to dTheta/dt and multiply rho
730  cell_src(i,j,k,RhoTheta_comp) += cell_data(i,j,k,Rho_comp) * dTdt * iexner;
731  }
732  });
733  }
734  } // mfi
735  } // OMP
736 }
void ApplySpongeZoneBCsForCC(const SpongeChoice &spongeChoice, const Geometry geom, const Box &bx, const Array4< Real > &cell_rhs, const Array4< const Real > &cell_data, const Array4< const Real > &r0, const Array4< const Real > &z_phys_cc)
Definition: ERF_ApplySpongeZoneBCs.cpp:7
void ApplySurfaceTreatment_BulkCoeff_CC(const Box &bx, const Array4< Real > &cell_rhs, const Array4< const Real > &cons_state, const Array4< const Real > &z_phys_cc, const Array4< const Real > &surface_state_arr)
Definition: ERF_ApplySurfaceTreatment_BulkCoeff.cpp:60
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:20
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
@ thetabar
Definition: ERF_DataStruct.H:97
@ m_y
Definition: ERF_DataStruct.H:25
@ m_x
Definition: ERF_DataStruct.H:24
DirectionSelector< 2 > ZDir
Definition: ERF_DirectionSelector.H:38
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenRTh(const amrex::Real rhotheta, const amrex::Real rdOcp, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:156
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:50
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
void NumericalDiffusion_Scal(const Box &bx, const int start_comp, const int num_comp, 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:18
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
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11
Definition: ERF_PlaneAverage.H:14
@ r0_comp
Definition: ERF_IndexDefines.H:63
@ cons
Definition: ERF_IndexDefines.H:158
@ theta
Definition: ERF_MM5.H:20
@ qv
Definition: ERF_Kessler.H:28
@ nc
Definition: ERF_Morrison.H:44
@ 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
real(c_double), private rhoi
Definition: ERF_module_mp_morr_two_moment.F90:188
bool rayleigh_damp_T
Definition: ERF_DampingStruct.H:87
amrex::Real rayleigh_dampcoef
Definition: ERF_DampingStruct.H:88
RayleighDampingType rayleigh_damping_type
Definition: ERF_DampingStruct.H:101
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > theta_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::Real if_init_surf_temp
Definition: ERF_DataStruct.H:1113
amrex::Real if_surf_temp_flux
Definition: ERF_DataStruct.H:1112
DampingChoice dampingChoice
Definition: ERF_DataStruct.H:1058
bool spatial_moisture_forcing
Definition: ERF_DataStruct.H:1145
RadiationType rad_type
Definition: ERF_DataStruct.H:1175
amrex::Real if_z0
Definition: ERF_DataStruct.H:1111
bool do_theta_advection
Definition: ERF_DataStruct.H:1140
bool custom_rhotheta_forcing
Definition: ERF_DataStruct.H:1137
bool custom_w_subsidence
Definition: ERF_DataStruct.H:1139
bool nudging_from_input_sounding
Definition: ERF_DataStruct.H:1149
amrex::Real if_Cd_scalar
Definition: ERF_DataStruct.H:1108
bool immersed_forcing_substep
Definition: ERF_DataStruct.H:1104
amrex::Real if_Olen_in
Definition: ERF_DataStruct.H:1115
bool use_num_diff
Definition: ERF_DataStruct.H:1164
bool four_stream_radiation
Definition: ERF_DataStruct.H:1101
bool custom_moisture_forcing
Definition: ERF_DataStruct.H:1138
amrex::Real num_diff_coeff
Definition: ERF_DataStruct.H:1165
amrex::Real if_surf_heating_rate
Definition: ERF_DataStruct.H:1114
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1060
MoistureType moisture_type
Definition: ERF_DataStruct.H:1171
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
PerturbationType pert_type
Definition: ERF_DataStruct.H:1161
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 spatial_rhotheta_forcing
Definition: ERF_DataStruct.H:1144
int ave_plane
Definition: ERF_DataStruct.H:1186
std::string sponge_type
Definition: ERF_SpongeStruct.H:58
Definition: ERF_TurbStruct.H:42
bool use_tke
Definition: ERF_TurbStruct.H:430
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:465
amrex::Vector< amrex::MultiFab > pb_cell
Definition: ERF_TurbPertStruct.H:640
void apply_tpi(const int &lev, const amrex::Box &vbx, const int &comp, const amrex::IndexType &m_ixtype, const amrex::Array4< amrex::Real > &src_arr, const amrex::Array4< amrex::Real const > &pert_cell)
Definition: ERF_TurbPertStruct.H:324
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: