112 BL_PROFILE_REGION(
"realbdy_compute_interior_ghost_RHS()");
120 Real F1 = 1./(nudge_factor*delta_t);
123 Real dT = bdy_time_interval;
125 int n_time =
static_cast<int>( (time-start_bdy_time) / dT);
126 int n_time_p1 = n_time + 1;
127 Real alpha = ((time-start_bdy_time) - n_time * dT) / dT;
130 if (time >= final_bdy_time) {
131 n_time =
static_cast<int>( (final_bdy_time - start_bdy_time)/ dT);
145 FArrayBox U_xlo, U_xhi, U_ylo, U_yhi;
146 FArrayBox V_xlo, V_xhi, V_ylo, V_yhi;
147 FArrayBox T_xlo, T_xhi, T_ylo, T_yhi;
168 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
169 int ivar_idx = var_map[ivar];
170 Box domain = geom.Domain();
171 auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
172 domain.convert(ixtype);
175 IntVect ng_vect(1,1,0);
176 Box gdom(domain); gdom.grow(ng_vect);
177 Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
185 U_xlo.resize(bx_xlo,1,The_Async_Arena()); U_xhi.resize(bx_xhi,1,The_Async_Arena());
186 U_ylo.resize(bx_ylo,1,The_Async_Arena()); U_yhi.resize(bx_yhi,1,The_Async_Arena());
187 }
else if (ivar == ivarV) {
188 V_xlo.resize(bx_xlo,1,The_Async_Arena()); V_xhi.resize(bx_xhi,1,The_Async_Arena());
189 V_ylo.resize(bx_ylo,1,The_Async_Arena()); V_yhi.resize(bx_yhi,1,The_Async_Arena());
190 }
else if (ivar == ivarT){
191 T_xlo.resize(bx_xlo,1,The_Async_Arena()); T_xhi.resize(bx_xhi,1,The_Async_Arena());
192 T_ylo.resize(bx_ylo,1,The_Async_Arena()); T_yhi.resize(bx_yhi,1,The_Async_Arena());
207 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
208 int ivar_idx = var_map[ivar];
209 Box domain = geom.Domain();
210 auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
211 domain.convert(ixtype);
212 const auto&
dom_lo = lbound(domain);
213 const auto&
dom_hi = ubound(domain);
216 #pragma omp parallel if (Gpu::notInLaunchRegion())
218 for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
220 IntVect ng_vect(1,1,0);
221 Box gtbx = grow(mfi.tilebox(ixtype.toIntVect()),ng_vect);
222 Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
228 Array4<Real> arr_xlo; Array4<Real> arr_xhi;
229 Array4<Real> arr_ylo; Array4<Real> arr_yhi;
231 arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
232 arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
233 }
else if (ivar == ivarV) {
234 arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
235 arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
236 }
else if (ivar == ivarT){
237 arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
238 arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
244 const auto& bdatxlo_n = bdy_data_xlo[n_time ][ivar].const_array();
245 const auto& bdatxlo_np1 = bdy_data_xlo[n_time_p1][ivar].const_array();
246 const auto& bdatxhi_n = bdy_data_xhi[n_time ][ivar].const_array();
247 const auto& bdatxhi_np1 = bdy_data_xhi[n_time_p1][ivar].const_array();
248 const auto& bdatylo_n = bdy_data_ylo[n_time ][ivar].const_array();
249 const auto& bdatylo_np1 = bdy_data_ylo[n_time_p1][ivar].const_array();
250 const auto& bdatyhi_n = bdy_data_yhi[n_time ][ivar].const_array();
251 const auto& bdatyhi_np1 = bdy_data_yhi[n_time_p1][ivar].const_array();
261 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
263 int ii = std::max(i ,
dom_lo.x);
265 int jj = std::max(j ,
dom_lo.y);
266 jj = std::min(jj,
dom_hi.y);
270 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
271 }
else if (ivar==ivarV) {
272 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
274 rho_interp = r_arr(i,j,k);
276 arr_xlo(i,j,k) = rho_interp * ( oma * bdatxlo_n (ii,jj,k,0)
277 +
alpha * bdatxlo_np1(ii,jj,k,0) );
279 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
282 ii = std::min(ii,
dom_hi.x);
283 int jj = std::max(j ,
dom_lo.y);
284 jj = std::min(jj,
dom_hi.y);
288 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
289 }
else if (ivar==ivarV) {
290 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
292 rho_interp = r_arr(i,j,k);
294 arr_xhi(i,j,k) = rho_interp * ( oma * bdatxhi_n (ii,jj,k,0)
295 +
alpha * bdatxhi_np1(ii,jj,k,0) );
299 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
301 int ii = std::max(i ,
dom_lo.x);
302 ii = std::min(ii,
dom_hi.x);
303 int jj = std::max(j ,
dom_lo.y);
308 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
309 }
else if (ivar==ivarV) {
310 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
312 rho_interp = r_arr(i,j,k);
314 arr_ylo(i,j,k) = rho_interp * ( oma * bdatylo_n (ii,jj,k,0)
315 +
alpha * bdatylo_np1(ii,jj,k,0) );
317 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
319 int ii = std::max(i ,
dom_lo.x);
320 ii = std::min(ii,
dom_hi.x);
322 jj = std::min(jj,
dom_hi.y);
326 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
327 }
else if (ivar==ivarV) {
328 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
330 rho_interp = r_arr(i,j,k);
332 arr_yhi(i,j,k) = rho_interp * ( oma * bdatyhi_n (ii,jj,k,0)
333 +
alpha * bdatyhi_np1(ii,jj,k,0) );
341 auto dx = geom.CellSizeArray();
342 auto ProbLo = geom.ProbLoArray();
343 auto ProbHi = geom.ProbHiArray();
344 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
345 int ivar_idx = ivar_map[ivar];
346 int icomp = comp_map[ivar];
348 Box domain = geom.Domain();
349 domain.convert(S_cur_data[ivar_idx].boxArray().ixType());
353 #pragma omp parallel if (Gpu::notInLaunchRegion())
355 for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
356 Box tbx = mfi.tilebox();
357 Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
363 Array4<Real> rhs_arr; Array4<Real> data_arr;
364 Array4<Real> arr_xlo; Array4<Real> arr_xhi;
365 Array4<Real> arr_ylo; Array4<Real> arr_yhi;
367 arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
368 arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
371 }
else if (ivar == ivarV) {
372 arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
373 arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
376 }
else if (ivar == ivarT){
377 arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
378 arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
385 Array4<Real> u_xlo = U_xlo.array(); Array4<Real> u_xhi = U_xhi.array();
386 Array4<Real> v_xlo = V_xlo.array(); Array4<Real> v_xhi = V_xhi.array();
387 Array4<Real> v_ylo = V_ylo.array(); Array4<Real> v_yhi = V_yhi.array();
390 width,
dx, ProbLo, ProbHi, F1, geom.Domain(),
391 tbx_xlo , tbx_xhi , tbx_ylo , tbx_yhi ,
392 arr_xlo , arr_xhi , arr_ylo , arr_yhi ,
393 u_xlo, u_xhi, v_xlo, v_xhi, v_ylo, v_yhi,
394 data_arr, rhs_arr, do_upwind);
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
const auto dx
Definition: ERF_InitCustomPertVels_ParticleTests.H:15
amrex::Real alpha
Definition: ERF_InitCustomPert_IsentropicVortex.H:8
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
void realbdy_interior_bxs_xy(const Box &bx, const Box &domain, const int &width, Box &bx_xlo, Box &bx_xhi, Box &bx_ylo, Box &bx_yhi, const IntVect &ng_vect, const bool get_int_ng)
Definition: ERF_InteriorGhostCells.cpp:23
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_relaxation(const int &icomp, const int &num_var, const int &width, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dx, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbLo, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbHi, const amrex::Real &F1, const amrex::Box &domain_cc, const amrex::Box &bx_xlo, const amrex::Box &bx_xhi, const amrex::Box &bx_ylo, const amrex::Box &bx_yhi, const amrex::Array4< const amrex::Real > &arr_xlo, const amrex::Array4< const amrex::Real > &arr_xhi, const amrex::Array4< const amrex::Real > &arr_ylo, const amrex::Array4< const amrex::Real > &arr_yhi, const amrex::Array4< const amrex::Real > &u_xlo, const amrex::Array4< const amrex::Real > &u_xhi, const amrex::Array4< const amrex::Real > &v_xlo, const amrex::Array4< const amrex::Real > &v_xhi, const amrex::Array4< const amrex::Real > &v_ylo, const amrex::Array4< const amrex::Real > &v_yhi, const amrex::Array4< const amrex::Real > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr, const bool &do_upwind)
Definition: ERF_Utils.H:180
@ U
Definition: ERF_IndexDefines.H:108
@ NumTypes
Definition: ERF_IndexDefines.H:112
@ T
Definition: ERF_IndexDefines.H:110
@ V
Definition: ERF_IndexDefines.H:109
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ yvel
Definition: ERF_IndexDefines.H:142