Function for computing the advective tendency for the update equations for all scalars other than rho and (rho theta) This routine has explicit expressions for all cases (terrain or not) when the horizontal and vertical spatial orders are <= 2, and calls more specialized functions when either (or both) spatial order(s) is greater than 2.
141 auto dxInv = cellSizeInv[0], dyInv = cellSizeInv[1], dzInv = cellSizeInv[2];
143 const Box
xbx = surroundingNodes(bx,0);
144 const Box
ybx = surroundingNodes(bx,1);
145 const Box
zbx = surroundingNodes(bx,2);
154 Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
156 if ( bx.smallEnd(0) == domain.smallEnd(0)) { bx_xlo = makeSlab( bx,0,domain.smallEnd(0));}
159 if ( bx.bigEnd(0) == domain.bigEnd(0)) { bx_xhi = makeSlab( bx,0,domain.bigEnd(0) );}
162 if ( bx.smallEnd(1) == domain.smallEnd(1)) { bx_ylo = makeSlab( bx,1,domain.smallEnd(1));}
165 if ( bx.bigEnd(1) == domain.bigEnd(1)) { bx_yhi = makeSlab( bx,1,domain.bigEnd(1) );}
168 for (
int n(0); n<ncomp; ++n) {
169 const int cons_index = icomp + n;
177 ParallelFor(
xbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
179 const int prim_index = cons_index - 1;
180 const Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i-1,j,k,prim_index));
181 (flx_arr[0])(i,j,k) = avg_xmom(i,j,k) * prim_on_face;
183 ParallelFor(
ybx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
185 const int prim_index = cons_index - 1;
186 const Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i,j-1,k,prim_index));
187 (flx_arr[1])(i,j,k) = avg_ymom(i,j,k) * prim_on_face;
189 ParallelFor(
zbx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
191 const int prim_index = cons_index - 1;
192 const Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i,j,k-1,prim_index));
193 (flx_arr[2])(i,j,k) = avg_zmom(i,j,k) * prim_on_face;
198 switch(horiz_adv_type) {
200 AdvectionSrcForScalarsVert<CENTERED2>(bx, cons_index, flx_arr, cell_prim,
201 avg_xmom, avg_ymom, avg_zmom,
202 horiz_upw_frac, vert_upw_frac, vert_adv_type);
205 AdvectionSrcForScalarsVert<UPWIND3>(bx, cons_index, flx_arr, cell_prim,
206 avg_xmom, avg_ymom, avg_zmom,
207 horiz_upw_frac, vert_upw_frac, vert_adv_type);
210 AdvectionSrcForScalarsVert<UPWIND3SL>(bx, cons_index, flx_arr, cell_prim,
211 avg_xmom, avg_ymom, avg_zmom,
212 horiz_upw_frac, vert_upw_frac, vert_adv_type);
215 AdvectionSrcForScalarsVert<CENTERED4>(bx, cons_index, flx_arr, cell_prim,
216 avg_xmom, avg_ymom, avg_zmom,
217 horiz_upw_frac, vert_upw_frac, vert_adv_type);
220 AdvectionSrcForScalarsVert<UPWIND5>(bx, cons_index, flx_arr, cell_prim,
221 avg_xmom, avg_ymom, avg_zmom,
222 horiz_upw_frac, vert_upw_frac, vert_adv_type);
225 AdvectionSrcForScalarsVert<CENTERED6>(bx, cons_index, flx_arr, cell_prim,
226 avg_xmom, avg_ymom, avg_zmom,
227 horiz_upw_frac, vert_upw_frac, vert_adv_type);
230 AdvectionSrcForScalarsWrapper<WENO3,WENO3>(bx, cons_index, flx_arr, cell_prim,
231 avg_xmom, avg_ymom, avg_zmom,
232 horiz_upw_frac, vert_upw_frac);
235 AdvectionSrcForScalarsWrapper<WENO5,WENO5>(bx, cons_index, flx_arr, cell_prim,
236 avg_xmom, avg_ymom, avg_zmom,
237 horiz_upw_frac, vert_upw_frac);
240 AdvectionSrcForScalarsWrapper<WENO7,WENO7>(bx, cons_index, flx_arr, cell_prim,
241 avg_xmom, avg_ymom, avg_zmom,
242 horiz_upw_frac, vert_upw_frac);
245 AdvectionSrcForScalarsWrapper<WENO_Z3,WENO_Z3>(bx, cons_index, flx_arr, cell_prim,
246 avg_xmom, avg_ymom, avg_zmom,
247 horiz_upw_frac, vert_upw_frac);
250 AdvectionSrcForScalarsWrapper<WENO_MZQ3,WENO_MZQ3>(bx, cons_index, flx_arr, cell_prim,
251 avg_xmom, avg_ymom, avg_zmom,
252 horiz_upw_frac, vert_upw_frac);
255 AdvectionSrcForScalarsWrapper<WENO_Z5,WENO_Z5>(bx, cons_index, flx_arr, cell_prim,
256 avg_xmom, avg_ymom, avg_zmom,
257 horiz_upw_frac, vert_upw_frac);
260 AdvectionSrcForScalarsWrapper<WENO_Z7,WENO_Z7>(bx, cons_index, flx_arr, cell_prim,
261 avg_xmom, avg_ymom, avg_zmom,
262 horiz_upw_frac, vert_upw_frac);
265 AMREX_ASSERT_WITH_MESSAGE(
false,
"Unknown advection scheme!");
269 ParallelFor(bx, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
271 if (detJ(i,j,k) > 0.)
273 Real invdetJ = 1.0 / detJ(i,j,k);
274 Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
276 advectionSrc(i,j,k,cons_index) = - invdetJ * mfsq * (
277 ( (flx_arr[0])(i+1,j,k) - (flx_arr[0])(i,j,k) ) * dxInv +
278 ( (flx_arr[1])(i,j+1,k) - (flx_arr[1])(i,j,k) ) * dyInv +
279 ( (flx_arr[2])(i,j,k+1) - (flx_arr[2])(i,j,k) ) * dzInv );
281 advectionSrc(i,j,k) = 0.;
289 avg_xmom, avg_ymom, avg_zmom,
290 detJ, cellSizeInv, do_lo);
294 avg_xmom, avg_ymom, avg_zmom,
300 avg_xmom, avg_ymom, avg_zmom,
301 detJ, cellSizeInv, do_lo);
305 avg_xmom, avg_ymom, avg_zmom,
void AdvectionSrcForScalars(const Box &bx, const int icomp, const int ncomp, const Array4< const Real > &avg_xmom, const Array4< const Real > &avg_ymom, const Array4< const Real > &avg_zmom, const Array4< const Real > &cell_prim, const Array4< Real > &advectionSrc, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mf_mx, const Array4< const Real > &mf_my, const AdvType horiz_adv_type, const AdvType vert_adv_type, const Real horiz_upw_frac, const Real vert_upw_frac, const GpuArray< const Array4< Real >, AMREX_SPACEDIM > &flx_arr, const Box &domain, const BCRec *bc_ptr_h)
Definition: ERF_AdvectionSrcForState.cpp:120
void AdvectionSrcForOpenBC_Tangent_Cons(const amrex::Box &bx, const int &dir, const int &icomp, const int &ncomp, const amrex::Array4< amrex::Real > &cell_rhs, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &avg_xmom, const amrex::Array4< const amrex::Real > &avg_ymom, const amrex::Array4< const amrex::Real > &avg_zmom, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const bool do_lo=false)
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ open
Definition: ERF_IndexDefines.H:215