ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_DiffusionSrcForMom_EB.cpp File Reference
#include <AMReX.H>
#include <AMReX_EB_Slopes_K.H>
#include <ERF_EB.H>
#include <ERF_Diffusion.H>
#include <ERF_IndexDefines.H>
#include <ERF_EBSlopes.H>
#include <ERF_DiffStruct.H>
#include <ERF_EBStruct.H>
Include dependency graph for ERF_DiffusionSrcForMom_EB.cpp:

Functions

void DiffusionSrcForMom_EB (const MFIter &mfi, [[maybe_unused]] const Box &domain, const Box &bxx, const Box &bxy, const Box &bxz, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< Real > &rho_w_rhs, const Array4< const Real > &u_arr, const Array4< const Real > &v_arr, const Array4< const Real > &w_arr, const Array4< const Real > &tau11, const Array4< const Real > &tau22, const Array4< const Real > &tau33, const Array4< const Real > &tau12, const Array4< const Real > &tau13, const Array4< const Real > &tau23, const Array4< const Real > &tau_eb13, const Array4< const Real > &tau_eb23, const Real *dx_arr, const GpuArray< Real, AMREX_SPACEDIM > &dxInv, const Array4< const Real > &mf_mx, const Array4< const Real > &mf_ux, const Array4< const Real > &mf_vx, const Array4< const Real > &mf_my, const Array4< const Real > &mf_uy, const Array4< const Real > &mf_vy, const SolverChoice &solverChoice, const eb_ &ebfact, [[maybe_unused]] const BCRec *d_bcrec_ptr)
 

Function Documentation

◆ DiffusionSrcForMom_EB()

void DiffusionSrcForMom_EB ( const MFIter &  mfi,
[[maybe_unused] ] const Box &  domain,
const Box &  bxx,
const Box &  bxy,
const Box &  bxz,
const Array4< Real > &  rho_u_rhs,
const Array4< Real > &  rho_v_rhs,
const Array4< Real > &  rho_w_rhs,
const Array4< const Real > &  u_arr,
const Array4< const Real > &  v_arr,
const Array4< const Real > &  w_arr,
const Array4< const Real > &  tau11,
const Array4< const Real > &  tau22,
const Array4< const Real > &  tau33,
const Array4< const Real > &  tau12,
const Array4< const Real > &  tau13,
const Array4< const Real > &  tau23,
const Array4< const Real > &  tau_eb13,
const Array4< const Real > &  tau_eb23,
const Real dx_arr,
const GpuArray< Real, AMREX_SPACEDIM > &  dxInv,
const Array4< const Real > &  mf_mx,
const Array4< const Real > &  mf_ux,
const Array4< const Real > &  mf_vx,
const Array4< const Real > &  mf_my,
const Array4< const Real > &  mf_uy,
const Array4< const Real > &  mf_vy,
const SolverChoice solverChoice,
const eb_ ebfact,
[[maybe_unused] ] const BCRec *  d_bcrec_ptr 
)

Function for computing the momentum RHS for diffusion operator without terrain.

Parameters
[in]mfiMultiFab Iterator
[in]domaincomputational domain
[in]bxxnodal x box for x-mom
[in]bxynodal y box for y-mom
[in]bxznodal z box for z-mom
[out]rho_u_rhsRHS for x-mom
[out]rho_v_rhsRHS for y-mom
[out]rho_w_rhsRHS for z-mom
[in]tau1111 stress
[in]tau2222 stress
[in]tau3333 stress
[in]tau1212 stress
[in]tau1313 stress
[in]tau2323 stress
[in]dxInvinverse cell size array
[in]mf_mmap factor at cell center
[in]ebfactEB factories for cell- and face-centered variables
62 {
63  BL_PROFILE_VAR("DiffusionSrcForMom_EB()",DiffusionSrcForMom_EB);
64 
65  DiffChoice dc = solverChoice.diffChoice;
66  const bool l_use_constAlpha = ( dc.molec_diff_type == MolecDiffType::ConstantAlpha );
67  Real mu_eff = (l_use_constAlpha) ? two * dc.dynamic_viscosity / dc.rho0_trans
68  : two * dc.dynamic_viscosity;
69 
70  auto dxinv = dxInv[0], dyinv = dxInv[1], dzinv = dxInv[2];
71  Real dx = dx_arr[0], dy = dx_arr[1], dz = dx_arr[2];
72  Real vol = dx * dy * dz;
73 
74  EBChoice ebChoice = solverChoice.ebChoice;
75  const bool l_no_slip = (ebChoice.eb_boundary_type == EBBoundaryType::NoSlipWall);
76  const bool l_surface_layer = (ebChoice.eb_boundary_type == EBBoundaryType::SurfaceLayer);
77  const bool l_constraint_x = solverChoice.diffChoice.eb_diff_constraint_x;
78  const bool l_constraint_y = solverChoice.diffChoice.eb_diff_constraint_y;
79  const bool l_constraint_z = solverChoice.diffChoice.eb_diff_constraint_z;
80 
81  // int n = 0;
82 
83  // EB u-factory
84  Array4<const EBCellFlag> u_cellflg = (ebfact.get_u_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
85  Array4<const Real > u_volfrac = (ebfact.get_u_const_factory())->getVolFrac().const_array(mfi);
86  Array4<const Real > u_volcent = (ebfact.get_u_const_factory())->getCentroid().const_array(mfi);
87  Array4<const Real > u_afrac_x = (ebfact.get_u_const_factory())->getAreaFrac()[0]->const_array(mfi);
88  Array4<const Real > u_afrac_y = (ebfact.get_u_const_factory())->getAreaFrac()[1]->const_array(mfi);
89  Array4<const Real > u_afrac_z = (ebfact.get_u_const_factory())->getAreaFrac()[2]->const_array(mfi);
90  Array4<const Real > u_bcent = (ebfact.get_u_const_factory())->getBndryCent().const_array(mfi);
91  Array4<const Real > u_bnorm = (ebfact.get_u_const_factory())->getBndryNorm().const_array(mfi);
92 
93  // EB v-factory
94  Array4<const EBCellFlag> v_cellflg = (ebfact.get_v_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
95  Array4<const Real > v_volfrac = (ebfact.get_v_const_factory())->getVolFrac().const_array(mfi);
96  Array4<const Real > v_volcent = (ebfact.get_v_const_factory())->getCentroid().const_array(mfi);
97  Array4<const Real > v_afrac_x = (ebfact.get_v_const_factory())->getAreaFrac()[0]->const_array(mfi);
98  Array4<const Real > v_afrac_y = (ebfact.get_v_const_factory())->getAreaFrac()[1]->const_array(mfi);
99  Array4<const Real > v_afrac_z = (ebfact.get_v_const_factory())->getAreaFrac()[2]->const_array(mfi);
100  Array4<const Real > v_bcent = (ebfact.get_v_const_factory())->getBndryCent().const_array(mfi);
101  Array4<const Real > v_bnorm = (ebfact.get_v_const_factory())->getBndryNorm().const_array(mfi);
102 
103  // EB w-factory
104  Array4<const EBCellFlag> w_cellflg = (ebfact.get_w_const_factory())->getMultiEBCellFlagFab()[mfi].const_array();
105  Array4<const Real > w_volfrac = (ebfact.get_w_const_factory())->getVolFrac().const_array(mfi);
106  Array4<const Real > w_volcent = (ebfact.get_w_const_factory())->getCentroid().const_array(mfi);
107  Array4<const Real > w_afrac_x = (ebfact.get_w_const_factory())->getAreaFrac()[0]->const_array(mfi);
108  Array4<const Real > w_afrac_y = (ebfact.get_w_const_factory())->getAreaFrac()[1]->const_array(mfi);
109  Array4<const Real > w_afrac_z = (ebfact.get_w_const_factory())->getAreaFrac()[2]->const_array(mfi);
110  Array4<const Real > w_bcent = (ebfact.get_w_const_factory())->getBndryCent().const_array(mfi);
111  Array4<const Real > w_bnorm = (ebfact.get_w_const_factory())->getBndryNorm().const_array(mfi);
112 
113  // x-momentum
114  ParallelFor(bxx,
115  [=] AMREX_GPU_DEVICE (int i, int j, int k)
116  {
117  if (u_volfrac(i,j,k)>zero) {
118 
119  // Inv Jacobian
120  Real mfsq = mf_ux(i,j,0) * mf_uy(i,j,0);
121 
122  Real diffContrib = ( (tau11(i , j , k ) * u_afrac_x(i+1,j ,k )
123  - tau11(i-1, j , k ) * u_afrac_x(i ,j ,k ) ) * dxinv * mfsq
124  + (tau12(i , j+1, k ) * u_afrac_y(i ,j+1,k )
125  - tau12(i , j , k ) * u_afrac_y(i ,j ,k ) ) * dyinv * mfsq
126  + (tau13(i , j , k+1) * u_afrac_z(i ,j ,k+1)
127  - tau13(i , j , k ) * u_afrac_z(i ,j ,k )) * dzinv );
128  diffContrib /= u_volfrac(i,j,k);
129 
130  rho_u_rhs(i,j,k) -= diffContrib;
131 
132  if (!l_constraint_x && u_cellflg(i,j,k).isSingleValued()) {
133 
134  Real axm = u_afrac_x(i ,j ,k );
135  Real axp = u_afrac_x(i+1,j ,k );
136  Real aym = u_afrac_y(i ,j ,k );
137  Real ayp = u_afrac_y(i ,j+1,k );
138  Real azm = u_afrac_z(i ,j ,k );
139  Real azp = u_afrac_z(i ,j ,k+1);
140 
141  Real adx = (axm-axp) * dy * dz;
142  Real ady = (aym-ayp) * dx * dz;
143  Real adz = (azm-azp) * dx * dy;
144 
145  Real barea = std::sqrt(adx*adx + ady*ady + adz*adz);
146 
147  Real dudn = zero;
148 
149  if (l_no_slip) {
150 
151  const RealVect bcent_eb {u_bcent(i,j,k,0), u_bcent(i,j,k,1), u_bcent(i,j,k,2)};
152  const Real Dirichlet_u {zero};
153  const Real Dirichlet_v {zero};
154  const Real Dirichlet_w {zero};
155 
156  GpuArray<Real,AMREX_SPACEDIM> slopes_u;
157  GpuArray<Real,AMREX_SPACEDIM> slopes_v;
158  GpuArray<Real,AMREX_SPACEDIM> slopes_w;
159 
160  slopes_u = erf_calc_slopes_eb_Dirichlet ( dx, dy, dz, i, j, k, bcent_eb, Dirichlet_u, u_arr, u_volcent, u_cellflg);
161  slopes_v = erf_calc_slopes_eb_Dirichlet_staggered( Vars::xvel, Vars::yvel, dx, dy, dz, i, j, k, bcent_eb, Dirichlet_v, v_arr, v_volcent, v_cellflg);
162  slopes_w = erf_calc_slopes_eb_Dirichlet_staggered( Vars::xvel, Vars::zvel, dx, dy, dz, i, j, k, bcent_eb, Dirichlet_w, w_arr, w_volcent, w_cellflg);
163 
164  Real dudx = slopes_u[0];
165  Real dudy = slopes_u[1];
166  Real dudz = slopes_u[2];
167  Real dvdx = slopes_v[0];
168  Real dvdy = slopes_v[1];
169  Real dwdx = slopes_w[0];
170  Real dwdz = slopes_w[2];
171 
172  Real tau11_eb = ( dudx - ( dudx + dvdy + dwdz ) / three );
173  Real tau12_eb = myhalf * (dudy + dvdx);
174  Real tau13_eb = myhalf * (dudz + dwdx);
175 
176  dudn = - mu_eff * (u_bnorm(i,j,k,0) * tau11_eb + u_bnorm(i,j,k,1) * tau12_eb + u_bnorm(i,j,k,2) * tau13_eb);
177 
178  } else if (l_surface_layer) {
179 
180  dudn = - tau_eb13(i,j,k);
181  }
182 
183  rho_u_rhs(i,j,k) -= barea * dudn / (vol * u_volfrac(i,j,k));
184  }
185  }
186 
187  });
188 
189  // y-momentum
190  ParallelFor(bxy,
191  [=] AMREX_GPU_DEVICE (int i, int j, int k)
192  {
193  if (v_volfrac(i,j,k)>zero) {
194 
195  // Inv Jacobian
196  Real mfsq = mf_vx(i,j,0) * mf_vy(i,j,0);
197 
198  Real diffContrib = ( (tau12(i+1, j , k ) * v_afrac_x(i+1,j ,k )
199  - tau12(i , j , k ) * v_afrac_x(i ,j ,k ) ) * dxinv * mfsq
200  + (tau22(i , j , k ) * v_afrac_y(i ,j+1,k )
201  - tau22(i , j-1, k ) * v_afrac_y(i ,j ,k ) ) * dyinv * mfsq
202  + (tau23(i , j , k+1) * v_afrac_z(i ,j ,k+1)
203  - tau23(i , j , k ) * v_afrac_z(i ,j ,k ) ) * dzinv );
204  diffContrib /= v_volfrac(i,j,k);
205 
206  rho_v_rhs(i,j,k) -= diffContrib;
207 
208  if (!l_constraint_y && l_no_slip && v_cellflg(i,j,k).isSingleValued()) {
209 
210  Real axm = v_afrac_x(i ,j ,k );
211  Real axp = v_afrac_x(i+1,j ,k );
212  Real aym = v_afrac_y(i ,j ,k );
213  Real ayp = v_afrac_y(i ,j+1,k );
214  Real azm = v_afrac_z(i ,j ,k );
215  Real azp = v_afrac_z(i ,j ,k+1);
216 
217  Real adx = (axm-axp) * dy * dz;
218  Real ady = (aym-ayp) * dx * dz;
219  Real adz = (azm-azp) * dx * dy;
220 
221  Real barea = std::sqrt(adx*adx + ady*ady + adz*adz);
222 
223  Real dvdn = 0.0;
224 
225  if (l_no_slip) {
226 
227  const RealVect bcent_eb {v_bcent(i,j,k,0), v_bcent(i,j,k,1), v_bcent(i,j,k,2)};
228  const Real Dirichlet_u {zero};
229  const Real Dirichlet_v {zero};
230  const Real Dirichlet_w {zero};
231 
232  GpuArray<Real,AMREX_SPACEDIM> slopes_u;
233  GpuArray<Real,AMREX_SPACEDIM> slopes_v;
234  GpuArray<Real,AMREX_SPACEDIM> slopes_w;
235 
236  slopes_u = erf_calc_slopes_eb_Dirichlet_staggered( Vars::yvel, Vars::xvel, dx, dy, dz, i, j, k, bcent_eb, Dirichlet_u, u_arr, u_volcent, u_cellflg);
237  slopes_v = erf_calc_slopes_eb_Dirichlet ( dx, dy, dz, i, j, k, bcent_eb, Dirichlet_v, v_arr, v_volcent, v_cellflg);
238  slopes_w = erf_calc_slopes_eb_Dirichlet_staggered( Vars::yvel, Vars::zvel, dx, dy, dz, i, j, k, bcent_eb, Dirichlet_w, w_arr, w_volcent, w_cellflg);
239 
240  Real dudx = slopes_u[0];
241  Real dudy = slopes_u[1];
242  Real dvdx = slopes_v[0];
243  Real dvdy = slopes_v[1];
244  Real dvdz = slopes_v[2];
245  Real dwdy = slopes_w[1];
246  Real dwdz = slopes_w[2];
247 
248  Real tau22_eb = ( dvdy - ( dudx + dvdy + dwdz ) / three );
249  Real tau12_eb = myhalf * (dudy + dvdx);
250  Real tau23_eb = myhalf * (dvdz + dwdy);
251 
252  dvdn = - mu_eff * (v_bnorm(i,j,k,0) * tau12_eb + v_bnorm(i,j,k,1) * tau22_eb + v_bnorm(i,j,k,2) * tau23_eb);
253 
254  } else if (l_surface_layer) {
255 
256  dvdn = - tau_eb23(i,j,k);
257  }
258 
259  rho_v_rhs(i,j,k) -= barea * dvdn / (vol * v_volfrac(i,j,k));
260  }
261  }
262  });
263 
264  // z-momentum
265  ParallelFor(bxz,
266  [=] AMREX_GPU_DEVICE (int i, int j, int k)
267  {
268  if (w_volfrac(i,j,k)>zero) {
269 
270  // Inv Jacobian
271  Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
272 
273  Real diffContrib = ( (tau13(i+1, j , k ) * w_afrac_x(i+1,j ,k )
274  - tau13(i , j , k ) * w_afrac_x(i ,j ,k ) ) * dxinv * mfsq
275  + (tau23(i , j+1, k ) * w_afrac_y(i ,j+1,k )
276  - tau23(i , j , k ) * w_afrac_y(i ,j ,k ) ) * dyinv * mfsq
277  + (tau33(i , j , k ) * w_afrac_z(i ,j ,k+1)
278  - tau33(i , j , k-1) * w_afrac_z(i ,j ,k ) ) * dzinv );
279  diffContrib /= w_volfrac(i,j,k);
280  rho_w_rhs(i,j,k) -= diffContrib;
281 
282  if (!l_constraint_z && l_no_slip && w_cellflg(i,j,k).isSingleValued()) {
283 
284  if (l_no_slip) {
285 
286  Real axm = w_afrac_x(i ,j ,k );
287  Real axp = w_afrac_x(i+1,j ,k );
288  Real aym = w_afrac_y(i ,j ,k );
289  Real ayp = w_afrac_y(i ,j+1,k );
290  Real azm = w_afrac_z(i ,j ,k );
291  Real azp = w_afrac_z(i ,j ,k+1);
292 
293  Real adx = (axm-axp) * dy * dz;
294  Real ady = (aym-ayp) * dx * dz;
295  Real adz = (azm-azp) * dx * dy;
296 
297  Real barea = std::sqrt(adx*adx + ady*ady + adz*adz);
298 
299  const RealVect bcent_eb {w_bcent(i,j,k,0), w_bcent(i,j,k,1), w_bcent(i,j,k,2)};
300  const Real Dirichlet_u {zero};
301  const Real Dirichlet_v {zero};
302  const Real Dirichlet_w {zero};
303 
304  GpuArray<Real,AMREX_SPACEDIM> slopes_u;
305  GpuArray<Real,AMREX_SPACEDIM> slopes_v;
306  GpuArray<Real,AMREX_SPACEDIM> slopes_w;
307 
308  slopes_u = erf_calc_slopes_eb_Dirichlet_staggered( Vars::zvel, Vars::xvel, dx, dy, dz, i, j, k, bcent_eb, Dirichlet_u, u_arr, u_volcent, u_cellflg);
309  slopes_v = erf_calc_slopes_eb_Dirichlet_staggered( Vars::zvel, Vars::yvel, dx, dy, dz, i, j, k, bcent_eb, Dirichlet_v, v_arr, v_volcent, v_cellflg);
310  slopes_w = erf_calc_slopes_eb_Dirichlet ( dx, dy, dz, i, j, k, bcent_eb, Dirichlet_w, w_arr, w_volcent, w_cellflg);
311 
312  Real dudx = slopes_u[0];
313  Real dudz = slopes_u[2];
314  Real dvdy = slopes_v[1];
315  Real dvdz = slopes_v[2];
316  Real dwdx = slopes_w[0];
317  Real dwdy = slopes_w[1];
318  Real dwdz = slopes_w[2];
319 
320  Real tau33_eb = ( dwdz - ( dudx + dvdy + dwdz ) / three );
321  Real tau13_eb = myhalf * (dudz + dwdx);
322  Real tau23_eb = myhalf * (dvdz + dwdy);
323 
324  Real dwdn = -(w_bnorm(i,j,k,0) * tau13_eb + w_bnorm(i,j,k,1) * tau23_eb + w_bnorm(i,j,k,2) * tau33_eb);
325 
326  rho_w_rhs(i,j,k) -= mu_eff * barea * dwdn / (vol * w_volfrac(i,j,k));
327 
328  }
329  }
330  }
331  });
332 
333 }
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
void DiffusionSrcForMom_EB(const MFIter &mfi, [[maybe_unused]] const Box &domain, const Box &bxx, const Box &bxy, const Box &bxz, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< Real > &rho_w_rhs, const Array4< const Real > &u_arr, const Array4< const Real > &v_arr, const Array4< const Real > &w_arr, const Array4< const Real > &tau11, const Array4< const Real > &tau22, const Array4< const Real > &tau33, const Array4< const Real > &tau12, const Array4< const Real > &tau13, const Array4< const Real > &tau23, const Array4< const Real > &tau_eb13, const Array4< const Real > &tau_eb23, const Real *dx_arr, const GpuArray< Real, AMREX_SPACEDIM > &dxInv, const Array4< const Real > &mf_mx, const Array4< const Real > &mf_ux, const Array4< const Real > &mf_vx, const Array4< const Real > &mf_my, const Array4< const Real > &mf_uy, const Array4< const Real > &mf_vy, const SolverChoice &solverChoice, const eb_ &ebfact, [[maybe_unused]] const BCRec *d_bcrec_ptr)
Definition: ERF_DiffusionSrcForMom_EB.cpp:34
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > erf_calc_slopes_eb_Dirichlet(amrex::Real dx, amrex::Real dy, amrex::Real dz, int i, int j, int k, amrex::RealVect const &bcent_eb, amrex::Real const state_eb, amrex::Array4< amrex::Real const > const &state, amrex::Array4< amrex::Real const > const &ccent, amrex::Array4< amrex::EBCellFlag const > const &flag)
Definition: ERF_EBSlopes.H:20
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > erf_calc_slopes_eb_Dirichlet_staggered(int igrid_query, int igrid_data, amrex::Real dx, amrex::Real dy, amrex::Real dz, int i, int j, int k, amrex::RealVect const &bcent_eb, amrex::Real const state_eb, amrex::Array4< amrex::Real const > const &state, amrex::Array4< amrex::Real const > const &ccent, amrex::Array4< amrex::EBCellFlag const > const &flag)
Definition: ERF_EBSlopes.H:146
@ tau_eb23
Definition: ERF_EBStruct.H:16
@ tau_eb13
Definition: ERF_EBStruct.H:16
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
const Real dy
Definition: ERF_InitCustomPert_ABL.H:24
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
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);} } })
amrex::Real Real
Definition: ERF_ShocInterface.H:19
eb_aux_ const * get_w_const_factory() const noexcept
Definition: ERF_EB.H:52
eb_aux_ const * get_v_const_factory() const noexcept
Definition: ERF_EB.H:51
eb_aux_ const * get_u_const_factory() const noexcept
Definition: ERF_EB.H:50
@ xvel
Definition: ERF_IndexDefines.H:159
@ zvel
Definition: ERF_IndexDefines.H:161
@ yvel
Definition: ERF_IndexDefines.H:160
Definition: ERF_DiffStruct.H:19
amrex::Real rho0_trans
Definition: ERF_DiffStruct.H:91
bool eb_diff_constraint_z
Definition: ERF_DiffStruct.H:99
bool eb_diff_constraint_x
Definition: ERF_DiffStruct.H:97
bool eb_diff_constraint_y
Definition: ERF_DiffStruct.H:98
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:84
amrex::Real dynamic_viscosity
Definition: ERF_DiffStruct.H:96
Definition: ERF_EBStruct.H:23
DiffChoice diffChoice
Definition: ERF_DataStruct.H:1088
EBChoice ebChoice
Definition: ERF_DataStruct.H:1092
Here is the call graph for this function: