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