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_PlaneAverage.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
59 {
60  BL_PROFILE_REGION("erf_make_sources()");
61 
62  // *****************************************************************************
63  // Initialize source to zero since we re-compute it every RK stage
64  // *****************************************************************************
65  if (is_slow_step) {
66  source.setVal(0.);
67  } else {
68  source.setVal(0.0,Rho_comp,2);
69  }
70 
71  const bool l_use_ndiff = solverChoice.use_num_diff;
72 
73  TurbChoice tc = solverChoice.turbChoice[level];
74  const bool l_use_KE = tc.use_tke;
75  const bool l_diff_KE = tc.diffuse_tke_3D;
76 
77  const Box& domain = geom.Domain();
78 
79  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
80  const GpuArray<Real, AMREX_SPACEDIM> dx = geom.CellSizeArray();
81 
82  MultiFab r_hse (base_state, make_alias, BaseState::r0_comp , 1);
83  MultiFab th_hse (base_state, make_alias, BaseState::th0_comp , 1);
84  MultiFab qv_hse (base_state, make_alias, BaseState::qv0_comp , 1);
85 
86  Real* thetabar = d_rayleigh_ptrs_at_lev[Rayleigh::thetabar];
87 
88  // flags to apply certain source terms in substep call only
89  bool use_Rayleigh_fast = ( (solverChoice.dampingChoice.rayleigh_damping_type == RayleighDampingType::FastExplicit) ||
91  bool use_ImmersedForcing_fast = solverChoice.immersed_forcing_substep;
92 
93  // flag for a moisture model
94  bool has_moisture = (solverChoice.moisture_type != MoistureType::None);
95 
96  // *****************************************************************************
97  // Planar averages for subsidence terms
98  // *****************************************************************************
99  Table1D<Real> dptr_r_plane, dptr_t_plane, dptr_qv_plane, dptr_qc_plane;
100  TableData<Real, 1> r_plane_tab, t_plane_tab, qv_plane_tab, qc_plane_tab;
101  bool compute_averages = ( is_slow_step && (dptr_wbar_sub || solverChoice.nudging_from_input_sounding) );
102 
103  if (compute_averages)
104  {
105  // The plane averaging operates at fixed z not fixed height so is not correct for variable dz
106  AMREX_ALWAYS_ASSERT(solverChoice.mesh_type != MeshType::VariableDz);
107 
108  //
109  // The call to "compute_averages" currently does all the components in one call
110  // We can then extract each component separately with the "line_average" call
111  //
112  // We need just one ghost cell in the vertical
113  //
114  IntVect ng_c(S_data[IntVars::cons].nGrowVect()); ng_c[2] = 1;
115  //
116  // With no moisture we only (rho) and (rho theta); with moisture we also do qv and qc
117  // We use the alias here to control ncomp inside the PlaneAverage
118  //
119  int ncomp = (!has_moisture) ? 2 : RhoQ2_comp+1;
120  MultiFab cons(S_data[IntVars::cons], make_alias, 0, ncomp);
121 
122  PlaneAverage cons_ave(&cons, geom, solverChoice.ave_plane, ng_c);
123  cons_ave.compute_averages(ZDir(), cons_ave.field());
124 
125  int ncell = cons_ave.ncell_line();
126 
127  Gpu::HostVector< Real> r_plane_h(ncell);
128  Gpu::DeviceVector< Real> r_plane_d(ncell);
129 
130  Gpu::HostVector< Real> t_plane_h(ncell);
131  Gpu::DeviceVector< Real> t_plane_d(ncell);
132 
133  cons_ave.line_average(Rho_comp , r_plane_h);
134  cons_ave.line_average(RhoTheta_comp, t_plane_h);
135 
136  Gpu::copyAsync(Gpu::hostToDevice, r_plane_h.begin(), r_plane_h.end(), r_plane_d.begin());
137  Gpu::copyAsync(Gpu::hostToDevice, t_plane_h.begin(), t_plane_h.end(), t_plane_d.begin());
138 
139  Real* dptr_r = r_plane_d.data();
140  Real* dptr_t = t_plane_d.data();
141 
142  Box tdomain = domain; tdomain.grow(2,ng_c[2]);
143  r_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
144  t_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
145 
146  int offset = ng_c[2];
147 
148  dptr_r_plane = r_plane_tab.table();
149  dptr_t_plane = t_plane_tab.table();
150  ParallelFor(ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
151  {
152  dptr_r_plane(k-offset) = dptr_r[k];
153  dptr_t_plane(k-offset) = dptr_t[k];
154  });
155 
156  if (has_moisture)
157  {
158  Gpu::HostVector< Real> qv_plane_h(ncell), qc_plane_h(ncell);
159  Gpu::DeviceVector<Real> qv_plane_d(ncell), qc_plane_d(ncell);
160 
161  // Water vapor
162  cons_ave.line_average(RhoQ1_comp, qv_plane_h);
163  Gpu::copyAsync(Gpu::hostToDevice, qv_plane_h.begin(), qv_plane_h.end(), qv_plane_d.begin());
164 
165  // Cloud water
166  cons_ave.line_average(RhoQ2_comp, qc_plane_h);
167  Gpu::copyAsync(Gpu::hostToDevice, qc_plane_h.begin(), qc_plane_h.end(), qc_plane_d.begin());
168 
169  Real* dptr_qv = qv_plane_d.data();
170  Real* dptr_qc = qc_plane_d.data();
171 
172  qv_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
173  qc_plane_tab.resize({tdomain.smallEnd(2)}, {tdomain.bigEnd(2)});
174 
175  dptr_qv_plane = qv_plane_tab.table();
176  dptr_qc_plane = qc_plane_tab.table();
177  ParallelFor(ncell, [=] AMREX_GPU_DEVICE (int k) noexcept
178  {
179  dptr_qv_plane(k-offset) = dptr_qv[k];
180  dptr_qc_plane(k-offset) = dptr_qc[k];
181  });
182  }
183  }
184 
185  // *****************************************************************************
186  // Radiation flux vector for four stream approximation
187  // *****************************************************************************
188  // NOTE: The fluxes live on w-faces
189  int klo = domain.smallEnd(0);
190  int khi = domain.bigEnd(2);
191  int nk = khi - klo + 2;
192  Gpu::DeviceVector<Real> radiation_flux(nk,zero);
193  Gpu::DeviceVector<Real> q_integral(nk,zero);
194  Real* rad_flux = radiation_flux.data();
195  Real* q_int = q_integral.data();
196 
197  // *****************************************************************************
198  // Define source term for cell-centered conserved variables, from
199  // one user-defined source terms for (rho theta) and (rho q_t)
200  // two radiation for (rho theta)
201  // three Rayleigh damping for (rho theta)
202  // Real(4.) custom forcing for (rho theta) and (rho Q1)
203  // Real(5.) custom subsidence for (rho theta) and (rho Q1)
204  // Real(6.) numerical diffusion for (rho theta)
205  // Real(7.) sponging
206  // Real(8.) turbulent perturbation
207  // Real(9.) nudging towards input sounding values (only for theta)
208  // 10a. Immersed forcing for terrain
209  // 10b. Immersed forcing for buildings
210  // Real(11.) Four stream radiation source for (rho theta)
211  // *****************************************************************************
212 
213  // ***********************************************************************************************
214  // Add remaining source terms
215  // ***********************************************************************************************
216 #ifdef _OPENMP
217 #pragma omp parallel if (Gpu::notInLaunchRegion())
218 #endif
219  {
220  for ( MFIter mfi(S_data[IntVars::cons],TileNoZ()); mfi.isValid(); ++mfi)
221  {
222  Box bx = mfi.tilebox();
223 
224  const Array4<const Real>& cell_data = S_data[IntVars::cons].array(mfi);
225  const Array4<const Real>& cell_prim = S_prim.array(mfi);
226  const Array4<Real> & cell_src = source.array(mfi);
227 
228  const Array4<const Real>& r0 = r_hse.const_array(mfi);
229  const Array4<const Real>& th0 = th_hse.const_array(mfi);
230  const Array4<const Real>& qv0 = qv_hse.const_array(mfi);
231 
232  const Array4<const Real>& z_cc_arr = z_phys_cc->const_array(mfi);
233 
234  const Array4<const Real>& t_blank_arr = (terrain_blank) ? terrain_blank->const_array(mfi) :
235  Array4<const Real>{};
236 
237 
238  // *************************************************************************************
239  // two Add radiation source terms to (rho theta)
240  // *************************************************************************************
241  if (solverChoice.rad_type != RadiationType::None && is_slow_step) {
242  auto const& qheating_arr = qheating_rates->const_array(mfi);
243  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
244  {
245  // Short-wavelength and long-wavelength radiation source terms
246  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) );
247  });
248  }
249 
250 
251  // *************************************************************************************
252  // three Add Rayleigh damping for (rho theta)
253  // *************************************************************************************
254  Real dampcoef = solverChoice.dampingChoice.rayleigh_dampcoef;
255 
256  if (solverChoice.dampingChoice.rayleigh_damp_T) {
257  if ((is_slow_step && !use_Rayleigh_fast) || (!is_slow_step && use_Rayleigh_fast)) {
258  int n = RhoTheta_comp;
259  int nr = Rho_comp;
260  int np = PrimTheta_comp;
261 
262  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
263  {
264  Real theta = cell_prim(i,j,k,np);
265  Real sinesq = d_sinesq_at_lev[k];
266  cell_src(i, j, k, n) -= dampcoef*sinesq * (theta - thetabar[k]) * cell_data(i,j,k,nr);
267  });
268  }
269  }
270 
271  // *************************************************************************************
272  // Real(4.) Add custom forcing for (rho theta)
273  // *************************************************************************************
274  if (solverChoice.custom_rhotheta_forcing && is_slow_step) {
275  const int n = RhoTheta_comp;
276  auto const& rhotheta_src_arr = rhotheta_src->const_array(mfi);
277  if (solverChoice.spatial_rhotheta_forcing)
278  {
279  if (solverChoice.custom_forcing_prim_vars) {
280  const int nr = Rho_comp;
281  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
282  {
283  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhotheta_src_arr(i, j, k);
284  });
285  } else {
286  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
287  {
288  cell_src(i, j, k, n) += rhotheta_src_arr(i, j, k);
289  });
290  }
291  } else {
292  if (solverChoice.custom_forcing_prim_vars) {
293  const int nr = Rho_comp;
294  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
295  {
296  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhotheta_src_arr(0, 0, k);
297  });
298  } else {
299  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
300  {
301  cell_src(i, j, k, n) += rhotheta_src_arr(0, 0, k);
302  });
303  }
304  }
305  }
306 
307  // *************************************************************************************
308  // Real(4.) Add custom forcing for RhoQ1
309  // *************************************************************************************
310  if (solverChoice.custom_moisture_forcing && is_slow_step) {
311  const int n = RhoQ1_comp;
312  auto const& rhoqt_src_arr = rhoqt_src->const_array(mfi);
313  if (solverChoice.spatial_moisture_forcing)
314  {
315  if (solverChoice.custom_forcing_prim_vars) {
316  const int nr = Rho_comp;
317  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
318  {
319  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhoqt_src_arr(i, j, k);
320  });
321  } else {
322  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
323  {
324  cell_src(i, j, k, n) += rhoqt_src_arr(i, j, k);
325  });
326  }
327  } else {
328  if (solverChoice.custom_forcing_prim_vars) {
329  const int nr = Rho_comp;
330  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
331  {
332  cell_src(i, j, k, n) += cell_data(i,j,k,nr) * rhoqt_src_arr(0, 0, k);
333  });
334  } else {
335  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
336  {
337  cell_src(i, j, k, n) += rhoqt_src_arr(0, 0, k);
338  });
339  }
340  }
341  }
342 
343  // *************************************************************************************
344  // Real(5.) Add custom subsidence for (rho theta)
345  // *************************************************************************************
346  if (solverChoice.custom_w_subsidence && is_slow_step && solverChoice.do_theta_advection) {
347  const int n = RhoTheta_comp;
348  if (solverChoice.custom_forcing_prim_vars) {
349  const int nr = Rho_comp;
350  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
351  {
352  Real dzInv = (z_cc_arr) ? one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : myhalf*dxInv[2];
353  Real T_hi = dptr_t_plane(k+1) / dptr_r_plane(k+1);
354  Real T_lo = dptr_t_plane(k-1) / dptr_r_plane(k-1);
355  Real wbar_cc = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
356  cell_src(i, j, k, n) -= cell_data(i,j,k,nr) * wbar_cc * (T_hi - T_lo) * dzInv;
357  });
358  } else {
359  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
360  {
361  Real dzInv = (z_cc_arr) ? one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : myhalf*dxInv[2];
362  Real T_hi = dptr_t_plane(k+1) / dptr_r_plane(k+1);
363  Real T_lo = dptr_t_plane(k-1) / dptr_r_plane(k-1);
364  Real wbar_cc = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
365  cell_src(i, j, k, n) -= wbar_cc * (T_hi - T_lo) * dzInv;
366  });
367  }
368  }
369 
370  // *************************************************************************************
371  // Real(5.) Add custom subsidence for RhoQ1 and RhoQ2
372  // *************************************************************************************
373  if (solverChoice.custom_w_subsidence && (solverChoice.moisture_type != MoistureType::None) && is_slow_step) {
374  const int nv = RhoQ1_comp;
375  if (solverChoice.custom_forcing_prim_vars) {
376  const int nr = Rho_comp;
377  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
378  {
379  Real dzInv = (z_cc_arr) ? one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : myhalf*dxInv[2];
380  Real Qv_hi = dptr_qv_plane(k+1) / dptr_r_plane(k+1);
381  Real Qv_lo = dptr_qv_plane(k-1) / dptr_r_plane(k-1);
382  Real Qc_hi = dptr_qc_plane(k+1) / dptr_r_plane(k+1);
383  Real Qc_lo = dptr_qc_plane(k-1) / dptr_r_plane(k-1);
384  Real wbar_cc = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
385  cell_src(i, j, k, nv ) -= cell_data(i,j,k,nr) * wbar_cc * (Qv_hi - Qv_lo) * dzInv;
386  cell_src(i, j, k, nv+1) -= cell_data(i,j,k,nr) * wbar_cc * (Qc_hi - Qc_lo) * dzInv;
387  });
388  } else {
389  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
390  {
391  Real dzInv = (z_cc_arr) ? one/ (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1)) : myhalf*dxInv[2];
392  Real Qv_hi = dptr_qv_plane(k+1) / dptr_r_plane(k+1);
393  Real Qv_lo = dptr_qv_plane(k-1) / dptr_r_plane(k-1);
394  Real Qc_hi = dptr_qc_plane(k+1) / dptr_r_plane(k+1);
395  Real Qc_lo = dptr_qc_plane(k-1) / dptr_r_plane(k-1);
396  Real wbar_cc = myhalf * (dptr_wbar_sub[k] + dptr_wbar_sub[k+1]);
397  cell_src(i, j, k, nv ) -= wbar_cc * (Qv_hi - Qv_lo) * dzInv;
398  cell_src(i, j, k, nv+1) -= wbar_cc * (Qc_hi - Qc_lo) * dzInv;
399  });
400  }
401  }
402 
403  // *************************************************************************************
404  // Real(6.) Add numerical diffusion for rho and (rho theta)
405  // *************************************************************************************
406  if (l_use_ndiff && is_slow_step)
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, 0, 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, 1, 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, RhoKE_comp, 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 == SpongeType::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.use_source_perturbation(level) && 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 
520  // Note this has been converted to K / s when it was read in;
521  const Real surf_heating_rate = solverChoice.if_surf_heating_rate;
522 
523  const Real Olen_in = solverChoice.if_Olen_in;
524 
525  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
526  {
527  const Real t_blank = t_blank_arr(i, j, k);
528  const Real t_blank_above = t_blank_arr(i, j, k+1);
529  const Real ux_cc_2r = myhalf * (u(i ,j ,k+1) + u(i+1,j ,k+1));
530  const Real uy_cc_2r = myhalf * (v(i ,j ,k+1) + v(i ,j+1,k+1));
531  const Real h_windspeed2r = std::sqrt(ux_cc_2r * ux_cc_2r + uy_cc_2r * uy_cc_2r);
532 
533  const Real theta = cell_data(i,j,k ,RhoTheta_comp) / cell_data(i,j,k ,Rho_comp);
534  const Real theta_neighbor = cell_data(i,j,k+1,RhoTheta_comp) / cell_data(i,j,k+1,Rho_comp);
535 
536  // SURFACE TEMP AND HEATING/COOLING RATE
537  if (init_surf_temp > zero) {
538  if (t_blank > 0 && (t_blank_above == zero)) { // force to MOST value
539  const Real surf_temp = init_surf_temp + surf_heating_rate*time;
540  const Real bc_forcing_rt_srf = -(cell_data(i,j,k-1,Rho_comp) * surf_temp - cell_data(i,j,k-1,RhoTheta_comp));
541  cell_src(i, j, k-1, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf; // k-1
542  }
543  }
544 
545  // SURFACE HEAT FLUX
546  if (tflux != Real(1e-8)){
547  if (t_blank > 0 && (t_blank_above == zero)) { // force to MOST value
548  Real psi_m = zero;
549  Real psi_h = zero;
550  Real psi_h_neighbor = zero;
551  Real ustar = h_windspeed2r * kappa / (std::log((Real(1.5)) * dx_z / z0) - psi_m);
552  const Real Olen = -ustar * ustar * ustar * theta / (kappa * ggg * tflux + tiny);
553  const Real zeta = (myhalf) * dx_z / Olen;
554  const Real zeta_neighbor = (Real(1.5)) * dx_z / Olen;
555 
556  // similarity functions
557  psi_m = sfuns.calc_psi_m(zeta);
558  psi_h = sfuns.calc_psi_h(zeta);
559  psi_h_neighbor = sfuns.calc_psi_h(zeta_neighbor);
560  ustar = h_windspeed2r * kappa / (std::log((Real(1.5)) * dx_z / z0) - psi_m);
561 
562  // prevent some unphysical math
563  if (!(ustar > zero && !std::isnan(ustar))) { ustar = zero; }
564  if (!(ustar < two && !std::isnan(ustar))) { ustar = two; }
565  if (psi_h_neighbor > std::log(Real(1.5) * dx_z / z0)) { psi_h_neighbor = std::log(Real(1.5) * dx_z / z0); }
566  if (psi_h > std::log(myhalf * dx_z / z0)) { psi_h = std::log(myhalf * dx_z / z0); }
567 
568  // We do not know the actual temperature so use cell above
569  const Real thetastar = theta * ustar * ustar / (kappa * ggg * Olen);
570  const Real surf_temp = theta_neighbor - thetastar / kappa * (std::log((Real(1.5)) * dx_z / z0) - psi_h_neighbor);
571  const Real tTarget = surf_temp + thetastar / kappa * (std::log((myhalf) * dx_z / z0) - psi_h);
572 
573  const Real bc_forcing_rt = -(cell_data(i,j,k,Rho_comp) * tTarget - cell_data(i,j,k,RhoTheta_comp));
574  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt;
575  }
576  }
577 
578  // OBUKHOV LENGTH
579  if (Olen_in != Real(1e-8)){
580  if (t_blank > 0 && (t_blank_above == zero)) { // force to MOST value
581  const Real Olen = Olen_in;
582  const Real zeta = (myhalf) * dx_z / Olen;
583  const Real zeta_neighbor = (Real(1.5)) * dx_z / Olen;
584 
585  // similarity functions
586  const Real psi_m = sfuns.calc_psi_m(zeta);
587  const Real psi_h = sfuns.calc_psi_h(zeta);
588  const Real psi_h_neighbor = sfuns.calc_psi_h(zeta_neighbor);
589  const Real ustar = h_windspeed2r * kappa / (std::log((Real(1.5)) * dx_z / z0) - psi_m);
590 
591  // We do not know the actual temperature so use cell above
592  const Real thetastar = theta * ustar * ustar / (kappa * ggg * Olen);
593  const Real surf_temp = theta_neighbor - thetastar / kappa * (std::log((Real(1.5)) * dx_z / z0) - psi_h_neighbor);
594  const Real tTarget = surf_temp + thetastar / kappa * (std::log((myhalf) * dx_z / z0) - psi_h);
595 
596  const Real bc_forcing_rt = -(cell_data(i,j,k,Rho_comp) * tTarget - cell_data(i,j,k,RhoTheta_comp));
597  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt;
598  }
599  }
600 
601  });
602  }
603 
604  // *************************************************************************************
605  // 10b. Add immersed source terms for buildings
606  // *************************************************************************************
607  if ((solverChoice.buildings_type == BuildingsType::ImmersedForcing ) &&
608  ((is_slow_step && !use_ImmersedForcing_fast) || (!is_slow_step && use_ImmersedForcing_fast)))
609  {
610  // geometric properties
611  const Real* dx_arr = geom.CellSize();
612  const Real dx_x = dx_arr[0];
613  const Real dx_y = dx_arr[1];
614 
615  const Real alpha_h = solverChoice.if_Cd_scalar;
616  const Real U_s = one; // unit velocity scale
617  const Real min_t_blank = Real(0.005);
618 
619  const Real init_surf_temp = solverChoice.if_init_surf_temp;
620 
621  // Note this has been converted to K / s when it was read in;
622  const Real surf_heating_rate = solverChoice.if_surf_heating_rate;
623 
624  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
625  {
626  Real t_blank = t_blank_arr(i, j, k);
627  Real t_blank_below = t_blank_arr(i, j, k-1);
628  Real t_blank_above = t_blank_arr(i, j, k+1);
629  Real t_blank_north = t_blank_arr(i , j+1, k);
630  Real t_blank_south = t_blank_arr(i , j-1, k);
631  Real t_blank_east = t_blank_arr(i+1, j , k);
632  Real t_blank_west = t_blank_arr(i-1, j , k);
633  if (t_blank < min_t_blank) { t_blank = zero; } // deal with situations where very small volfrac exist
634  if (t_blank_below < min_t_blank) { t_blank_below = zero; }
635  if (t_blank_north < min_t_blank) { t_blank_north = zero; }
636  if (t_blank_south < min_t_blank) { t_blank_south = zero; }
637  if (t_blank_east < min_t_blank) { t_blank_east = zero; }
638  if (t_blank_west < min_t_blank) { t_blank_west = zero; }
639 
640  Real dx_z = (z_cc_arr) ? (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-1)) : dx[2];
641  Real drag_coefficient = alpha_h / std::pow(dx_x*dx_y*dx_z, one/three);
642 
643  // SURFACE TEMP AND HEATING/COOLING RATE
644  if (init_surf_temp > zero) {
645  const Real surf_temp = init_surf_temp + surf_heating_rate*time;
646  if (t_blank > 0 && (t_blank_above == zero) && (t_blank_below == one)) { // building roof
647  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
648  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
649 
650  } else if (((t_blank > 0 && t_blank < t_blank_west && t_blank_east == zero) ||
651  (t_blank > 0 && t_blank < t_blank_east && t_blank_west == zero) ||
652  (t_blank > 0 && t_blank < t_blank_north && t_blank_south == zero) ||
653  (t_blank > 0 && t_blank < t_blank_south && t_blank_north == zero))) {
654  // this should enter for just building walls
655  // walls are currently separated to allow for flexible in the future to heat walls differently
656 
657  // south face
658  if ((t_blank < t_blank_north) && (t_blank_north == one)) {
659  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
660  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
661  }
662 
663  // north face
664  if ((t_blank < t_blank_south) && (t_blank_south == one)) {
665  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
666  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
667  }
668 
669  // west face
670  if ((t_blank < t_blank_east) && (t_blank_east == one)) {
671  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
672  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
673  }
674 
675  // east face
676  if ((t_blank < t_blank_west) && (t_blank_west == one)) {
677  const Real bc_forcing_rt_srf = -(cell_data(i,j,k,Rho_comp) * surf_temp - cell_data(i,j,k,RhoTheta_comp));
678  cell_src(i, j, k, RhoTheta_comp) -= drag_coefficient * U_s * bc_forcing_rt_srf;
679  }
680 
681  }
682  }
683  });
684  }
685 
686  // *************************************************************************************
687  // Real(11.) Add 4 stream radiation src to RhoTheta
688  // *************************************************************************************
689  if (solverChoice.four_stream_radiation && has_moisture && is_slow_step)
690  {
691  AMREX_ALWAYS_ASSERT((bx.smallEnd(2) == klo) && (bx.bigEnd(2) == khi));
692  Real D = Real(3.75e-6); // [s^-1]
693  Real F0 = 70; // [W/m^2]
694  Real F1 = 22; // [W/m^2]
695  Real krad = 85; // [m^2 kg^-1]
696  Real qt_i = Real(0.008);
697 
698  Box xybx = makeSlab(bx,2,klo);
699  ParallelFor(xybx, [=] AMREX_GPU_DEVICE(int i, int j, int /*k*/) noexcept
700  {
701  // Inclusive scan at w-faces for the Q integral (also find "i" values)
702  q_int[0] = zero;
703  Real zi = myhalf * (z_cc_arr(i,j,khi) + z_cc_arr(i,j,khi-1));
704  Real rhoi = myhalf * (cell_data(i,j,khi,Rho_comp) + cell_data(i,j,khi-1,Rho_comp));
705  for (int k(klo+1); k<=khi+1; ++k) {
706  int lk = k - klo;
707  // Average to w-faces when looping w-faces
708  Real dz = (z_cc_arr) ? myhalf * (z_cc_arr(i,j,k) - z_cc_arr(i,j,k-2)) : dx[2];
709  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;
710  Real qt_hi = cell_data(i,j,k ,RhoQ1_comp) + cell_data(i,j,k ,RhoQ2_comp);
711  Real qt_lo = cell_data(i,j,k-1,RhoQ1_comp) + cell_data(i,j,k-1,RhoQ2_comp);
712  if ( (qt_lo > qt_i) && (qt_hi < qt_i) ) {
713  zi = myhalf * (z_cc_arr(i,j,k) + z_cc_arr(i,j,k-1));
714  rhoi = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i,j,k-1,Rho_comp));
715  }
716  }
717 
718  // Decompose the integral to get the fluxes at w-faces
719  Real q_int_inf = q_int[khi+1];
720  for (int k(klo); k<=khi+1; ++k) {
721  int lk = k - klo;
722  Real z = myhalf * (z_cc_arr(i,j,k) + z_cc_arr(i,j,k-1));
723  rad_flux[lk] = F1*std::exp(-q_int[lk]) + F0*std::exp(-(q_int_inf - q_int[lk]));
724  if (z > zi) {
725  rad_flux[lk] += rhoi * Cp_d * D * ( std::pow(z-zi,Real(4.)/three)/Real(4.) + zi*std::pow(z-zi,one/three) ) ;
726  }
727  }
728 
729  // Compute the radiative heating source
730  for (int k(klo); k<=khi; ++k) {
731  int lk = k - klo;
732  // Average to w-faces when looping CC
733  Real dzInv = (z_cc_arr) ? one/ (myhalf * (z_cc_arr(i,j,k+1) - z_cc_arr(i,j,k-1))) : dxInv[2];
734  // NOTE: Fnet = Up - Dn (all fluxes are up here)
735  // dT/dt = dF/dz * (1/(-rho*Cp))
736  Real dTdt = (rad_flux[lk+1] - rad_flux[lk]) * dzInv / (-cell_data(i,j,k,Rho_comp)*Cp_d);
737  Real qv = cell_data(i,j,k,RhoQ1_comp)/cell_data(i,j,k,Rho_comp);
738  Real iexner = one/getExnergivenRTh(cell_data(i,j,k,RhoTheta_comp), R_d/Cp_d, qv);
739  // Convert dT/dt to dTheta/dt and multiply rho
740  cell_src(i,j,k,RhoTheta_comp) += cell_data(i,j,k,Rho_comp) * dTdt * iexner;
741  }
742  });
743  }
744  } // mfi
745  } // OMP
746 }
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:39
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:31
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:40
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
@ 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:55
#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)
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
void NumericalDiffusion_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:77
@ th0_comp
Definition: ERF_IndexDefines.H:76
@ r0_comp
Definition: ERF_IndexDefines.H:73
@ cons
Definition: ERF_IndexDefines.H:192
@ theta
Definition: ERF_MM5.H:20
@ qv
Definition: ERF_Kessler.H:29
@ nr
Definition: ERF_Morrison.H:45
@ xvel
Definition: ERF_IndexDefines.H:175
@ cons
Definition: ERF_IndexDefines.H:174
@ yvel
Definition: ERF_IndexDefines.H:176
@ dz
Definition: ERF_AdvanceWSM6.cpp:104
@ zi
Definition: ERF_AdvanceWSM6.cpp:133
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:1208
static MeshType mesh_type
Definition: ERF_DataStruct.H:1134
amrex::Real if_surf_temp_flux
Definition: ERF_DataStruct.H:1207
DampingChoice dampingChoice
Definition: ERF_DataStruct.H:1144
bool spatial_moisture_forcing
Definition: ERF_DataStruct.H:1243
RadiationType rad_type
Definition: ERF_DataStruct.H:1303
amrex::Real if_z0
Definition: ERF_DataStruct.H:1206
bool do_theta_advection
Definition: ERF_DataStruct.H:1238
bool custom_rhotheta_forcing
Definition: ERF_DataStruct.H:1235
bool custom_w_subsidence
Definition: ERF_DataStruct.H:1237
bool nudging_from_input_sounding
Definition: ERF_DataStruct.H:1247
amrex::Real if_Cd_scalar
Definition: ERF_DataStruct.H:1203
bool immersed_forcing_substep
Definition: ERF_DataStruct.H:1199
amrex::Real if_Olen_in
Definition: ERF_DataStruct.H:1210
bool use_num_diff
Definition: ERF_DataStruct.H:1295
bool four_stream_radiation
Definition: ERF_DataStruct.H:1196
bool custom_moisture_forcing
Definition: ERF_DataStruct.H:1236
amrex::Real num_diff_coeff
Definition: ERF_DataStruct.H:1296
amrex::Real if_surf_heating_rate
Definition: ERF_DataStruct.H:1209
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1146
MoistureType moisture_type
Definition: ERF_DataStruct.H:1299
bool custom_forcing_prim_vars
Definition: ERF_DataStruct.H:1241
static BuildingsType buildings_type
Definition: ERF_DataStruct.H:1128
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1125
SpongeChoice spongeChoice
Definition: ERF_DataStruct.H:1145
static InitType init_type
Definition: ERF_DataStruct.H:1119
bool use_source_perturbation(int lev) const
Definition: ERF_DataStruct.H:1274
bool hindcast_surface_bcs
Definition: ERF_DataStruct.H:1345
bool spatial_rhotheta_forcing
Definition: ERF_DataStruct.H:1242
int ave_plane
Definition: ERF_DataStruct.H:1313
static SpongeType sponge_type
Definition: ERF_SpongeStruct.H:81
Definition: ERF_TurbStruct.H:82
bool use_tke
Definition: ERF_TurbStruct.H:471
bool diffuse_tke_3D
Definition: ERF_TurbStruct.H:506
amrex::Vector< amrex::MultiFab > pb_cell
Definition: ERF_TurbPertStruct.H:667
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:351
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: