Function for computing the advective tendency for the momentum equations 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.
73 auto dxInv = cellSizeInv[0], dyInv = cellSizeInv[1], dzInv = cellSizeInv[2];
75 AMREX_ALWAYS_ASSERT(bxz.smallEnd(2) > 0);
78 Box box2d_u(bxx); box2d_u.setRange(2,0); box2d_u.grow({3,3,0});
79 Box box2d_v(bxy); box2d_v.setRange(2,0); box2d_v.grow({3,3,0});
80 FArrayBox mf_u_invFAB(box2d_u,1,The_Async_Arena());
81 FArrayBox mf_v_invFAB(box2d_v,1,The_Async_Arena());
82 const Array4<Real>& mf_u_inv = mf_u_invFAB.array();
83 const Array4<Real>& mf_v_inv = mf_v_invFAB.array();
85 ParallelFor(box2d_u, box2d_v,
86 [=] AMREX_GPU_DEVICE (
int i,
int j,
int) noexcept
88 mf_u_inv(i,j,0) = 1. / mf_u(i,j,0);
90 [=] AMREX_GPU_DEVICE (
int i,
int j,
int) noexcept
92 mf_v_inv(i,j,0) = 1. / mf_v(i,j,0);
96 amrex::ignore_unused(use_terrain);
97 ParallelFor(bxx, bxy, bxz,
98 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
100 rho_u_rhs(i, j, k) = 0.0;
102 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
104 rho_v_rhs(i, j, k) = 0.0;
106 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
108 rho_w_rhs(i, j, k) = 0.0;
115 ParallelFor(bxx, bxy, bxz,
116 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
118 Real xflux_hi = 0.25 * (rho_u(i, j , k) * mf_u_inv(i,j,0) + rho_u(i+1, j , k) * mf_u_inv(i+1,j,0)) * (u(i+1,j,k) + u(i,j,k));
119 Real xflux_lo = 0.25 * (rho_u(i, j , k) * mf_u_inv(i,j,0) + rho_u(i-1, j , k) * mf_u_inv(i-1,j,0)) * (u(i-1,j,k) + u(i,j,k));
121 Real yflux_hi = 0.25 * (rho_v(i, j+1, k) * mf_v_inv(i,j+1,0) + rho_v(i-1, j+1, k) * mf_v_inv(i-1,j+1,0)) * (u(i,j+1,k) + u(i,j,k));
122 Real yflux_lo = 0.25 * (rho_v(i, j , k) * mf_v_inv(i,j ,0) + rho_v(i-1, j , k) * mf_v_inv(i-1,j ,0)) * (u(i,j-1,k) + u(i,j,k));
124 Real zflux_hi = 0.25 * (Omega(i, j, k+1) + Omega(i-1, j, k+1)) * (u(i,j,k+1) + u(i,j,k));
125 Real zflux_lo = 0.25 * (Omega(i, j, k ) + Omega(i-1, j, k )) * (u(i,j,k-1) + u(i,j,k));
127 Real mfsq = mf_u(i,j,0) * mf_u(i,j,0);
129 Real advectionSrc = (xflux_hi - xflux_lo) * dxInv * mfsq
130 + (yflux_hi - yflux_lo) * dyInv * mfsq
131 + (zflux_hi - zflux_lo) * dzInv;
132 rho_u_rhs(i, j, k) = -advectionSrc;
134 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
136 Real xflux_hi = 0.25 * (rho_u(i+1, j, k) * mf_u_inv(i+1,j,0) + rho_u(i+1, j-1, k) * mf_u_inv(i+1,j-1,0)) * (v(i+1,j,k) + v(i,j,k));
137 Real xflux_lo = 0.25 * (rho_u(i , j, k) * mf_u_inv(i ,j,0) + rho_u(i , j-1, k) * mf_u_inv(i ,j-1,0)) * (v(i-1,j,k) + v(i,j,k));
139 Real yflux_hi = 0.25 * (rho_v(i ,j+1,k) * mf_v_inv(i,j+1,0) + rho_v(i ,j ,k) * mf_v_inv(i,j ,0)) * (v(i,j+1,k) + v(i,j,k));
140 Real yflux_lo = 0.25 * (rho_v(i ,j ,k) * mf_v_inv(i,j ,0) + rho_v(i ,j-1,k) * mf_v_inv(i,j-1,0) ) * (v(i,j-1,k) + v(i,j,k));
142 Real zflux_hi = 0.25 * (Omega(i, j, k+1) + Omega(i, j-1, k+1)) * (v(i,j,k+1) + v(i,j,k));
143 Real zflux_lo = 0.25 * (Omega(i, j, k ) + Omega(i, j-1, k )) * (v(i,j,k-1) + v(i,j,k));
145 Real mfsq = mf_v(i,j,0) * mf_v(i,j,0);
147 Real advectionSrc = (xflux_hi - xflux_lo) * dxInv * mfsq
148 + (yflux_hi - yflux_lo) * dyInv * mfsq
149 + (zflux_hi - zflux_lo) * dzInv;
150 rho_v_rhs(i, j, k) = -advectionSrc;
152 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
154 Real xflux_hi = 0.25*(rho_u(i+1,j ,k) + rho_u(i+1, j, k-1)) * mf_u_inv(i+1,j ,0) * (w(i+1,j,k) + w(i,j,k));
155 Real xflux_lo = 0.25*(rho_u(i ,j ,k) + rho_u(i , j, k-1)) * mf_u_inv(i ,j ,0) * (w(i-1,j,k) + w(i,j,k));
157 Real yflux_hi = 0.25*(rho_v(i ,j+1,k) + rho_v(i, j+1, k-1)) * mf_v_inv(i ,j+1,0) * (w(i,j+1,k) + w(i,j,k));
158 Real yflux_lo = 0.25*(rho_v(i ,j ,k) + rho_v(i, j , k-1)) * mf_v_inv(i ,j ,0) * (w(i,j-1,k) + w(i,j,k));
160 Real zflux_lo = 0.25 * (Omega(i,j,k) + Omega(i,j,k-1)) * (w(i,j,k) + w(i,j,k-1));
162 Real zflux_hi = (k == hi_z_face) ? Omega(i,j,k) * w(i,j,k) :
163 0.25 * (Omega(i,j,k) + Omega(i,j,k+1)) * (w(i,j,k) + w(i,j,k+1));
165 Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);
167 Real advectionSrc = (xflux_hi - xflux_lo) * dxInv * mfsq
168 + (yflux_hi - yflux_lo) * dyInv * mfsq
169 + (zflux_hi - zflux_lo) * dzInv;
170 rho_w_rhs(i, j, k) = -advectionSrc;
175 AdvectionSrcForMomVert_N<CENTERED2>(bxx, bxy, bxz,
176 rho_u_rhs, rho_v_rhs, rho_w_rhs,
177 rho_u, rho_v, Omega, u, v, w,
180 horiz_upw_frac, vert_upw_frac,
181 vert_adv_type, lo_z_face, hi_z_face);
183 AdvectionSrcForMomVert_N<UPWIND3>(bxx, bxy, bxz,
184 rho_u_rhs, rho_v_rhs, rho_w_rhs,
185 rho_u, rho_v, Omega, u, v, w,
188 horiz_upw_frac, vert_upw_frac,
189 vert_adv_type, lo_z_face, hi_z_face);
191 AdvectionSrcForMomVert_N<CENTERED4>(bxx, bxy, bxz,
192 rho_u_rhs, rho_v_rhs, rho_w_rhs,
193 rho_u, rho_v, Omega, u, v, w,
196 horiz_upw_frac, vert_upw_frac,
197 vert_adv_type, lo_z_face, hi_z_face);
199 AdvectionSrcForMomVert_N<UPWIND5>(bxx, bxy, bxz,
200 rho_u_rhs, rho_v_rhs, rho_w_rhs,
201 rho_u, rho_v, Omega, u, v, w,
204 horiz_upw_frac, vert_upw_frac,
205 vert_adv_type, lo_z_face, hi_z_face);
207 AdvectionSrcForMomVert_N<CENTERED6>(bxx, bxy, bxz,
208 rho_u_rhs, rho_v_rhs, rho_w_rhs,
209 rho_u, rho_v, Omega, u, v, w,
212 horiz_upw_frac, vert_upw_frac,
213 vert_adv_type, lo_z_face, hi_z_face);
215 AMREX_ASSERT_WITH_MESSAGE(
false,
"Unknown advection scheme!");
224 ParallelFor(bxx, bxy, bxz,
225 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
227 Real xflux_hi = 0.25 * (rho_u(i,j,k) * mf_u_inv(i,j,0) + rho_u(i+1,j,k) * mf_u_inv(i+1,j,0)) *
228 (u(i+1,j,k) + u(i,j,k)) * 0.5 * (ax(i,j,k) + ax(i+1,j,k));
230 Real xflux_lo = 0.25 * (rho_u(i,j,k) * mf_u_inv(i,j,0) + rho_u(i-1,j,k) * mf_u_inv(i-1,j,0)) *
231 (u(i-1,j,k) + u(i,j,k)) * 0.5 * (ax(i,j,k) + ax(i-1,j,k));
234 Real yflux_hi = 0.25 * (rho_v(i,j+1,k)*mf_v_inv(i,j+1,0) + rho_v(i-1,j+1,k)*mf_v_inv(i-1,j+1,0)) *
235 (u(i,j+1,k) + u(i,j,k)) * met_h_zeta_yhi;
238 Real yflux_lo = 0.25 * (rho_v(i,j ,k)*mf_v_inv(i,j ,0) + rho_v(i-1,j ,k)*mf_v_inv(i-1,j ,0)) *
239 (u(i,j-1,k) + u(i,j,k)) * met_h_zeta_ylo;
241 Real zflux_hi = 0.25 * (Omega(i,j,k+1) + Omega(i-1,j,k+1)) * (u(i,j,k+1) + u(i,j,k)) *
242 0.5 * (az(i,j,k+1) + az(i-1,j,k+1));
243 Real zflux_lo = 0.25 * (Omega(i,j,k ) + Omega(i-1,j,k )) * (u(i,j,k-1) + u(i,j,k)) *
244 0.5 * (az(i,j,k ) + az(i-1,j,k ));
246 Real mfsq = mf_u(i,j,0) * mf_u(i,j,0);
248 Real advectionSrc = (xflux_hi - xflux_lo) * dxInv * mfsq
249 + (yflux_hi - yflux_lo) * dyInv * mfsq
250 + (zflux_hi - zflux_lo) * dzInv;
252 rho_u_rhs(i, j, k) = -advectionSrc / (0.5 * (detJ(i,j,k) + detJ(i-1,j,k)));
254 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
258 Real xflux_hi = 0.25 * (rho_u(i+1,j,k)*mf_u_inv(i+1,j,0) + rho_u(i+1,j-1,k)*mf_u_inv(i+1,j-1,0)) *
259 (v(i+1,j,k) + v(i,j,k)) * met_h_zeta_xhi;
262 Real xflux_lo = 0.25 * (rho_u(i, j, k)*mf_u_inv(i ,j,0) + rho_u(i ,j-1,k)*mf_u_inv(i-1,j ,0)) *
263 (v(i-1,j,k) + v(i,j,k)) * met_h_zeta_xlo;
265 Real yflux_hi = 0.25 * (rho_v(i,j+1,k)*mf_v_inv(i,j+1,0) + rho_v(i,j ,k) * mf_v_inv(i,j ,0)) *
266 (v(i,j+1,k) + v(i,j,k)) * 0.5 * (ay(i,j,k) + ay(i,j+1,k));
268 Real yflux_lo = 0.25 * (rho_v(i,j ,k)*mf_v_inv(i,j ,0) + rho_v(i,j-1,k) * mf_v_inv(i,j-1,0)) *
269 (v(i,j-1,k) + v(i,j,k)) * 0.5 * (ay(i,j,k) + ay(i,j-1,k));
271 Real zflux_hi = 0.25 * (Omega(i,j,k+1) + Omega(i, j-1, k+1)) * (v(i,j,k+1) + v(i,j,k)) *
272 0.5 * (az(i,j,k+1) + az(i,j-1,k+1));
273 Real zflux_lo = 0.25 * (Omega(i,j,k ) + Omega(i, j-1, k )) * (v(i,j,k-1) + v(i,j,k)) *
274 0.5 * (az(i,j,k ) + az(i,j-1,k ));
276 Real mfsq = mf_v(i,j,0) * mf_v(i,j,0);
278 Real advectionSrc = (xflux_hi - xflux_lo) * dxInv * mfsq
279 + (yflux_hi - yflux_lo) * dyInv * mfsq
280 + (zflux_hi - zflux_lo) * dzInv;
282 rho_v_rhs(i, j, k) = -advectionSrc / (0.5 * (detJ(i,j,k) + detJ(i,j-1,k)));
284 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
287 Real xflux_hi = 0.25*(rho_u(i+1,j ,k) + rho_u(i+1,j,k-1)) * mf_u_inv(i+1,j,0) *
288 (w(i+1,j,k) + w(i,j,k)) * met_h_zeta_xhi;
291 Real xflux_lo = 0.25*(rho_u(i ,j ,k) + rho_u(i ,j,k-1)) * mf_u_inv(i ,j,0) *
292 (w(i-1,j,k) + w(i,j,k)) * met_h_zeta_xlo;
295 Real yflux_hi = 0.25*(rho_v(i,j+1,k) + rho_v(i,j+1,k-1)) * mf_v_inv(i,j+1,0) *
296 (w(i,j+1,k) + w(i,j,k)) * met_h_zeta_yhi;
299 Real yflux_lo = 0.25*(rho_v(i,j ,k) + rho_v(i,j ,k-1)) * mf_v_inv(i,j ,0) *
300 (w(i,j-1,k) + w(i,j,k)) * met_h_zeta_ylo;
302 Real zflux_lo = 0.25 * (Omega(i,j,k) + Omega(i,j,k-1)) * (w(i,j,k) + w(i,j,k-1));
304 Real zflux_hi = (k == hi_z_face) ? Omega(i,j,k) * w(i,j,k) * az(i,j,k):
305 0.25 * (Omega(i,j,k) + Omega(i,j,k+1)) * (w(i,j,k) + w(i,j,k+1)) *
306 0.5 * (az(i,j,k) + az(i,j,k+1));
308 Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);
310 Real advectionSrc = (xflux_hi - xflux_lo) * dxInv * mfsq
311 + (yflux_hi - yflux_lo) * dyInv * mfsq
312 + (zflux_hi - zflux_lo) * dzInv;
314 rho_w_rhs(i, j, k) = -advectionSrc / (0.5*(detJ(i,j,k) + detJ(i,j,k-1)));
319 AdvectionSrcForMomVert<CENTERED2>(bxx, bxy, bxz,
320 rho_u_rhs, rho_v_rhs, rho_w_rhs,
321 rho_u, rho_v, Omega, u, v, w, z_nd, ax, ay, az, detJ,
322 cellSizeInv, mf_m, mf_u_inv, mf_v_inv,
323 horiz_upw_frac, vert_upw_frac,
324 vert_adv_type, lo_z_face, hi_z_face);
326 AdvectionSrcForMomVert<UPWIND3>(bxx, bxy, bxz,
327 rho_u_rhs, rho_v_rhs, rho_w_rhs,
328 rho_u, rho_v, Omega, u, v, w, z_nd, ax, ay, az, detJ,
329 cellSizeInv, mf_m, mf_u_inv, mf_v_inv,
330 horiz_upw_frac, vert_upw_frac,
331 vert_adv_type, lo_z_face, hi_z_face);
333 AdvectionSrcForMomVert<CENTERED4>(bxx, bxy, bxz,
334 rho_u_rhs, rho_v_rhs, rho_w_rhs,
335 rho_u, rho_v, Omega, u, v, w, z_nd, ax, ay, az, detJ,
336 cellSizeInv, mf_m, mf_u_inv, mf_v_inv,
337 horiz_upw_frac, vert_upw_frac,
338 vert_adv_type, lo_z_face, hi_z_face);
340 AdvectionSrcForMomVert<UPWIND5>(bxx, bxy, bxz,
341 rho_u_rhs, rho_v_rhs, rho_w_rhs,
342 rho_u, rho_v, Omega, u, v, w, z_nd, ax, ay, az, detJ,
343 cellSizeInv, mf_m, mf_u_inv, mf_v_inv,
344 horiz_upw_frac, vert_upw_frac,
345 vert_adv_type, lo_z_face, hi_z_face);
347 AdvectionSrcForMomVert<CENTERED6>(bxx, bxy, bxz,
348 rho_u_rhs, rho_v_rhs, rho_w_rhs,
349 rho_u, rho_v, Omega, u, v, w, z_nd, ax, ay, az, detJ,
350 cellSizeInv, mf_m, mf_u_inv, mf_v_inv,
351 horiz_upw_frac, vert_upw_frac,
352 vert_adv_type, lo_z_face, hi_z_face);
354 AMREX_ASSERT_WITH_MESSAGE(
false,
"Unknown advection scheme!");
367 Box tbx(surroundingNodes(bx,0));
368 Box tby(surroundingNodes(bx,1));
369 Box tbz(surroundingNodes(bx,2)); tbz.growLo(2,-1); tbz.growHi(2,-1);
371 const int domhi_z = domain.bigEnd(2);
376 Box tbx_xlo, tby_xlo, tbz_xlo;
377 if (tbx.smallEnd(0) == domain.smallEnd(0)) { tbx_xlo = makeSlab(tbx,0,domain.smallEnd(0));}
378 if (tby.smallEnd(0) == domain.smallEnd(0)) { tby_xlo = makeSlab(tby,0,domain.smallEnd(0));}
379 if (tbz.smallEnd(0) == domain.smallEnd(0)) { tbz_xlo = makeSlab(tbz,0,domain.smallEnd(0));}
386 ay, az, detJ, cellSizeInv,
390 ax, ay, az, detJ, cellSizeInv,
395 Box tbx_xhi, tby_xhi, tbz_xhi;
396 if (tbx.bigEnd(0) == domain.bigEnd(0)+1) { tbx_xhi = makeSlab(tbx,0,domain.bigEnd(0)+1);}
397 if (tby.bigEnd(0) == domain.bigEnd(0)) { tby_xhi = makeSlab(tby,0,domain.bigEnd(0) );}
398 if (tbz.bigEnd(0) == domain.bigEnd(0)) { tbz_xhi = makeSlab(tbz,0,domain.bigEnd(0) );}
403 ay, az, detJ, cellSizeInv);
406 ax, ay, az, detJ, cellSizeInv,
411 Box tbx_ylo, tby_ylo, tbz_ylo;
412 if (tbx.smallEnd(1) == domain.smallEnd(1)) { tbx_ylo = makeSlab(tbx,1,domain.smallEnd(1));}
413 if (tby.smallEnd(1) == domain.smallEnd(1)) { tby_ylo = makeSlab(tby,1,domain.smallEnd(1));}
414 if (tbz.smallEnd(1) == domain.smallEnd(1)) { tbz_ylo = makeSlab(tbz,1,domain.smallEnd(1));}
419 ax, az, detJ, cellSizeInv,
424 ax, ay, az, detJ, cellSizeInv,
429 Box tbx_yhi, tby_yhi, tbz_yhi;
430 if (tbx.bigEnd(1) == domain.bigEnd(1)) { tbx_yhi = makeSlab(tbx,1,domain.bigEnd(1) );}
431 if (tby.bigEnd(1) == domain.bigEnd(1)+1) { tby_yhi = makeSlab(tby,1,domain.bigEnd(1)+1);}
432 if (tbz.bigEnd(1) == domain.bigEnd(1)) { tbz_yhi = makeSlab(tbz,1,domain.bigEnd(1) );}
436 ax, az, detJ, cellSizeInv);
440 ax, ay, az, detJ, cellSizeInv,
void AdvectionSrcForOpenBC_Tangent_Ymom(const amrex::Box &bxy, const int &dir, const amrex::Array4< amrex::Real > &rho_v_rhs, const amrex::Array4< const amrex::Real > &v, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &Omega, const amrex::Array4< const amrex::Real > &ay, const amrex::Array4< const amrex::Real > &az, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const bool do_lo=false)
void AdvectionSrcForOpenBC_Tangent_Xmom(const amrex::Box &bxx, const int &dir, const amrex::Array4< amrex::Real > &rho_u_rhs, const amrex::Array4< const amrex::Real > &u, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &Omega, const amrex::Array4< const amrex::Real > &ax, const amrex::Array4< const amrex::Real > &az, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const bool do_lo=false)
void AdvectionSrcForOpenBC_Normal(const amrex::Box &bx, const int &dir, const amrex::Array4< amrex::Real > &rhs_arr, const amrex::Array4< const amrex::Real > &vel_norm_arr, const amrex::Array4< const amrex::Real > &cell_data_arr, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const bool do_lo=false)
void AdvectionSrcForOpenBC_Tangent_Zmom(const amrex::Box &bxz, const int &dir, const amrex::Array4< amrex::Real > &rho_w_rhs, const amrex::Array4< const amrex::Real > &w, const amrex::Array4< const amrex::Real > &rho_u, const amrex::Array4< const amrex::Real > &rho_v, const amrex::Array4< const amrex::Real > &Omega, const amrex::Array4< const amrex::Real > &ax, const amrex::Array4< const amrex::Real > &ay, const amrex::Array4< const amrex::Real > &az, const amrex::Array4< const amrex::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const int domhi_z, const bool do_lo=false)
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterJ(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:266
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterK(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:221
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterI(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:310
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ open
Definition: ERF_IndexDefines.H:186