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.
145 auto dxInv = cellSizeInv[0], dyInv = cellSizeInv[1], dzInv = cellSizeInv[2];
147 const Box xbx = surroundingNodes(bx,0);
148 const Box ybx = surroundingNodes(bx,1);
149 const Box zbx = surroundingNodes(bx,2);
158 Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
160 if ( bx.smallEnd(0) == domain.smallEnd(0)) { bx_xlo = makeSlab( bx,0,domain.smallEnd(0));}
163 if ( bx.bigEnd(0) == domain.bigEnd(0)) { bx_xhi = makeSlab( bx,0,domain.bigEnd(0) );}
166 if ( bx.smallEnd(1) == domain.smallEnd(1)) { bx_ylo = makeSlab( bx,1,domain.smallEnd(1));}
169 if ( bx.bigEnd(1) == domain.bigEnd(1)) { bx_yhi = makeSlab( bx,1,domain.bigEnd(1) );}
178 ParallelFor(xbx, ncomp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
180 const int cons_index = icomp + n;
181 const int prim_index = cons_index - 1;
182 const Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i-1,j,k,prim_index));
183 (flx_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * prim_on_face;
185 ParallelFor(ybx, ncomp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
187 const int cons_index = icomp + n;
188 const int prim_index = cons_index - 1;
189 const Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i,j-1,k,prim_index));
190 (flx_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * prim_on_face;
192 ParallelFor(zbx, ncomp,[=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
194 const int cons_index = icomp + n;
195 const int prim_index = cons_index - 1;
196 const Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i,j,k-1,prim_index));
197 (flx_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * prim_on_face;
202 switch(horiz_adv_type) {
204 AdvectionSrcForScalarsVert<CENTERED2>(bx, ncomp, icomp, flx_arr, cell_prim,
205 avg_xmom, avg_ymom, avg_zmom,
206 horiz_upw_frac, vert_upw_frac, vert_adv_type);
209 AdvectionSrcForScalarsVert<UPWIND3>(bx, ncomp, icomp, flx_arr, cell_prim,
210 avg_xmom, avg_ymom, avg_zmom,
211 horiz_upw_frac, vert_upw_frac, vert_adv_type);
214 AdvectionSrcForScalarsVert<CENTERED4>(bx, ncomp, icomp, flx_arr, cell_prim,
215 avg_xmom, avg_ymom, avg_zmom,
216 horiz_upw_frac, vert_upw_frac, vert_adv_type);
219 AdvectionSrcForScalarsVert<UPWIND5>(bx, ncomp, icomp, flx_arr, cell_prim,
220 avg_xmom, avg_ymom, avg_zmom,
221 horiz_upw_frac, vert_upw_frac, vert_adv_type);
224 AdvectionSrcForScalarsVert<CENTERED6>(bx, ncomp, icomp, flx_arr, cell_prim,
225 avg_xmom, avg_ymom, avg_zmom,
226 horiz_upw_frac, vert_upw_frac, vert_adv_type);
229 AdvectionSrcForScalarsWrapper<WENO3,WENO3>(bx, ncomp, icomp, flx_arr, cell_prim,
230 avg_xmom, avg_ymom, avg_zmom,
231 horiz_upw_frac, vert_upw_frac);
234 AdvectionSrcForScalarsWrapper<WENO5,WENO5>(bx, ncomp, icomp, flx_arr, cell_prim,
235 avg_xmom, avg_ymom, avg_zmom,
236 horiz_upw_frac, vert_upw_frac);
239 AdvectionSrcForScalarsWrapper<WENO7,WENO7>(bx, ncomp, icomp, flx_arr, cell_prim,
240 avg_xmom, avg_ymom, avg_zmom,
241 horiz_upw_frac, vert_upw_frac);
244 AdvectionSrcForScalarsWrapper<WENO_Z3,WENO_Z3>(bx, ncomp, icomp, flx_arr, cell_prim,
245 avg_xmom, avg_ymom, avg_zmom,
246 horiz_upw_frac, vert_upw_frac);
249 AdvectionSrcForScalarsWrapper<WENO_MZQ3,WENO_MZQ3>(bx, ncomp, icomp, flx_arr, cell_prim,
250 avg_xmom, avg_ymom, avg_zmom,
251 horiz_upw_frac, vert_upw_frac);
254 AdvectionSrcForScalarsWrapper<WENO_Z5,WENO_Z5>(bx, ncomp, icomp, flx_arr, cell_prim,
255 avg_xmom, avg_ymom, avg_zmom,
256 horiz_upw_frac, vert_upw_frac);
259 AdvectionSrcForScalarsWrapper<WENO_Z7,WENO_Z7>(bx, ncomp, icomp, flx_arr, cell_prim,
260 avg_xmom, avg_ymom, avg_zmom,
261 horiz_upw_frac, vert_upw_frac);
264 AMREX_ASSERT_WITH_MESSAGE(
false,
"Unknown advection scheme!");
281 ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
283 const int cons_index = icomp + n;
284 (flx_tmp_arr[0])(i,j,k,cons_index) = (flx_arr[0])(i,j,k,cons_index);
285 (flx_tmp_arr[1])(i,j,k,cons_index) = (flx_arr[1])(i,j,k,cons_index);
286 (flx_tmp_arr[2])(i,j,k,cons_index) = (flx_arr[2])(i,j,k,cons_index);
290 ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
292 const int cons_index = icomp + n;
293 const int prim_index = cons_index - 1;
295 Real max_val = max_s_ptr[cons_index];
296 Real min_val = min_s_ptr[cons_index];
298 Real invdetJ = (detJ(i,j,k) > 0.) ? 1. / detJ(i,j,k) : 1.;
299 Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);
301 Real RHS = - invdetJ * mfsq * (
302 ( (flx_arr[0])(i+1,j ,k ,cons_index) - (flx_arr[0])(i,j,k,cons_index) ) * dxInv +
303 ( (flx_arr[1])(i ,j+1,k ,cons_index) - (flx_arr[1])(i,j,k,cons_index) ) * dyInv +
304 ( (flx_arr[2])(i ,j ,k+1,cons_index) - (flx_arr[2])(i,j,k,cons_index) ) * dzInv );
311 Real tmp_upd = cur_cons(i,j,k,cons_index) + RHS*dt;
312 if (tmp_upd<min_val || tmp_upd>max_val) {
314 if (avg_xmom(i+1,j,k)>0.0) {
315 (flx_tmp_arr[0])(i+1,j,k,cons_index) = avg_xmom(i+1,j,k) * cell_prim(i ,j,k,prim_index);
317 (flx_tmp_arr[0])(i+1,j,k,cons_index) = avg_xmom(i+1,j,k) * cell_prim(i+1,j,k,prim_index);
319 if (avg_ymom(i,j+1,k)>0.0) {
320 (flx_tmp_arr[1])(i,j+1,k,cons_index) = avg_ymom(i,j+1,k) * cell_prim(i,j ,k,prim_index);
322 (flx_tmp_arr[1])(i,j+1,k,cons_index) = avg_ymom(i,j+1,k) * cell_prim(i,j+1,k,prim_index);
324 if (avg_zmom(i,j,k+1)>0.0) {
325 (flx_tmp_arr[2])(i,j,k+1,cons_index) = avg_zmom(i,j,k+1) * cell_prim(i,j,k ,prim_index);
327 (flx_tmp_arr[2])(i,j,k+1,cons_index) = avg_zmom(i,j,k+1) * cell_prim(i,j,k+1,prim_index);
331 if (avg_xmom(i,j,k)>0.0) {
332 (flx_tmp_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * cell_prim(i-1,j,k,prim_index);
334 (flx_tmp_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * cell_prim(i ,j,k,prim_index);
336 if (avg_ymom(i,j,k)>0.0) {
337 (flx_tmp_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * cell_prim(i,j-1,k,prim_index);
339 (flx_tmp_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * cell_prim(i,j ,k,prim_index);
341 if (avg_zmom(i,j,k)>0.0) {
342 (flx_tmp_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * cell_prim(i,j,k-1,prim_index);
344 (flx_tmp_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * cell_prim(i,j,k ,prim_index);
350 ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
352 const int cons_index = icomp + n;
353 (flx_arr[0])(i,j,k,cons_index) = (flx_tmp_arr[0])(i,j,k,cons_index);
354 (flx_arr[1])(i,j,k,cons_index) = (flx_tmp_arr[1])(i,j,k,cons_index);
355 (flx_arr[2])(i,j,k,cons_index) = (flx_tmp_arr[2])(i,j,k,cons_index);
359 ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int n) noexcept
361 const int cons_index = icomp + n;
362 if (detJ(i,j,k) > 0.)
364 Real invdetJ = 1.0 / detJ(i,j,k);
365 Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);
367 advectionSrc(i,j,k,cons_index) = - invdetJ * mfsq * (
368 ( (flx_arr[0])(i+1,j,k,cons_index) - (flx_arr[0])(i ,j,k,cons_index) ) * dxInv +
369 ( (flx_arr[1])(i,j+1,k,cons_index) - (flx_arr[1])(i,j ,k,cons_index) ) * dyInv +
370 ( (flx_arr[2])(i,j,k+1,cons_index) - (flx_arr[2])(i,j,k ,cons_index) ) * dzInv );
372 advectionSrc(i,j,k,cons_index) = 0.;
380 avg_xmom, avg_ymom, avg_zmom,
381 detJ, cellSizeInv, do_lo);
385 avg_xmom, avg_ymom, avg_zmom,
391 avg_xmom, avg_ymom, avg_zmom,
392 detJ, cellSizeInv, do_lo);
396 avg_xmom, avg_ymom, avg_zmom,
void AdvectionSrcForScalars(const Real &dt, 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 > &cur_cons, const Array4< const Real > &cell_prim, const Array4< Real > &advectionSrc, const bool &use_mono_adv, Real *max_s_ptr, Real *min_s_ptr, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mf_m, 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 GpuArray< Array4< Real >, AMREX_SPACEDIM > &flx_tmp_arr, const Box &domain, const BCRec *bc_ptr_h)
Definition: ERF_AdvectionSrcForState.cpp:119
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:186