189 BL_PROFILE_REGION(
"realbdy_compute_interior_ghost_RHS()");
197 if (width > set_width+1) width -= 1;
200 Real F1 = 1./(10.*delta_t);
201 Real F2 = 1./(50.*delta_t);
204 Real dT = bdy_time_interval;
209 Real time_since_start = time;
211 int n_time =
static_cast<int>( time_since_start / dT);
212 Real alpha = (time_since_start - n_time * dT) / dT;
213 AMREX_ALWAYS_ASSERT( alpha >= 0. && alpha <= 1.0);
214 Real oma = 1.0 - alpha;
216 int n_time_p1 = n_time + 1;
217 if ((time == stop_time) && (alpha==0)) {
229 FArrayBox U_xlo, U_xhi, U_ylo, U_yhi;
230 FArrayBox V_xlo, V_xhi, V_ylo, V_yhi;
231 FArrayBox T_xlo, T_xhi, T_ylo, T_yhi;
252 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
253 int ivar_idx = var_map[ivar];
254 Box domain = geom.Domain();
255 auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
256 domain.convert(ixtype);
261 IntVect ng_vect{2,2,0};
262 Box gdom(domain); gdom.grow(ng_vect);
263 Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
271 U_xlo.resize(bx_xlo,1,The_Async_Arena()); U_xhi.resize(bx_xhi,1,The_Async_Arena());
272 U_ylo.resize(bx_ylo,1,The_Async_Arena()); U_yhi.resize(bx_yhi,1,The_Async_Arena());
273 }
else if (ivar == ivarV) {
274 V_xlo.resize(bx_xlo,1,The_Async_Arena()); V_xhi.resize(bx_xhi,1,The_Async_Arena());
275 V_ylo.resize(bx_ylo,1,The_Async_Arena()); V_yhi.resize(bx_yhi,1,The_Async_Arena());
276 }
else if (ivar == ivarT){
277 T_xlo.resize(bx_xlo,1,The_Async_Arena()); T_xhi.resize(bx_xhi,1,The_Async_Arena());
278 T_ylo.resize(bx_ylo,1,The_Async_Arena()); T_yhi.resize(bx_yhi,1,The_Async_Arena());
293 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
294 int ivar_idx = var_map[ivar];
295 Box domain = geom.Domain();
296 auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
297 domain.convert(ixtype);
298 const auto&
dom_lo = lbound(domain);
299 const auto&
dom_hi = ubound(domain);
302 #pragma omp parallel if (Gpu::notInLaunchRegion())
304 for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
307 IntVect ng_vect{2,2,0};
308 Box gtbx = grow(mfi.tilebox(ixtype.toIntVect()),ng_vect);
309 Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
315 Array4<Real> arr_xlo; Array4<Real> arr_xhi;
316 Array4<Real> arr_ylo; Array4<Real> arr_yhi;
318 arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
319 arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
320 }
else if (ivar == ivarV) {
321 arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
322 arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
323 }
else if (ivar == ivarT){
324 arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
325 arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
331 const auto& bdatxlo_n = bdy_data_xlo[n_time ][ivar].const_array();
332 const auto& bdatxlo_np1 = bdy_data_xlo[n_time_p1][ivar].const_array();
333 const auto& bdatxhi_n = bdy_data_xhi[n_time ][ivar].const_array();
334 const auto& bdatxhi_np1 = bdy_data_xhi[n_time_p1][ivar].const_array();
335 const auto& bdatylo_n = bdy_data_ylo[n_time ][ivar].const_array();
336 const auto& bdatylo_np1 = bdy_data_ylo[n_time_p1][ivar].const_array();
337 const auto& bdatyhi_n = bdy_data_yhi[n_time ][ivar].const_array();
338 const auto& bdatyhi_np1 = bdy_data_yhi[n_time_p1][ivar].const_array();
347 int offset = set_width - 1;
348 if (width > set_width)
offset = width;
351 ParallelFor(tbx_xlo, tbx_xhi,
352 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
354 int ii = std::max(i ,
dom_lo.x);
356 int jj = std::max(j ,
dom_lo.y);
357 jj = std::min(jj,
dom_hi.y);
361 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
362 }
else if (ivar==ivarV) {
363 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
365 rho_interp = r_arr(i,j,k);
367 arr_xlo(i,j,k) = rho_interp * ( oma * bdatxlo_n (ii,jj,k,0)
368 + alpha * bdatxlo_np1(ii,jj,k,0) );
370 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
373 ii = std::min(ii,
dom_hi.x);
374 int jj = std::max(j ,
dom_lo.y);
375 jj = std::min(jj,
dom_hi.y);
379 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
380 }
else if (ivar==ivarV) {
381 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
383 rho_interp = r_arr(i,j,k);
385 arr_xhi(i,j,k) = rho_interp * ( oma * bdatxhi_n (ii,jj,k,0)
386 + alpha * bdatxhi_np1(ii,jj,k,0) );
389 ParallelFor(tbx_ylo, tbx_yhi,
390 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
392 int ii = std::max(i ,
dom_lo.x);
393 ii = std::min(ii,
dom_hi.x);
394 int jj = std::max(j ,
dom_lo.y);
399 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
400 }
else if (ivar==ivarV) {
401 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
403 rho_interp = r_arr(i,j,k);
405 arr_ylo(i,j,k) = rho_interp * ( oma * bdatylo_n (ii,jj,k,0)
406 + alpha * bdatylo_np1(ii,jj,k,0) );
408 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept
410 int ii = std::max(i ,
dom_lo.x);
411 ii = std::min(ii,
dom_hi.x);
413 jj = std::min(jj,
dom_hi.y);
417 rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
418 }
else if (ivar==ivarV) {
419 rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
421 rho_interp = r_arr(i,j,k);
423 arr_yhi(i,j,k) = rho_interp * ( oma * bdatyhi_n (ii,jj,k,0)
424 + alpha * bdatyhi_np1(ii,jj,k,0) );
436 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
437 int ivar_idx = ivar_map[ivar];
438 int icomp = comp_map[ivar];
440 Box domain = geom.Domain();
441 auto ix_type = S_old_data[ivar_idx].boxArray().ixType();
442 auto iv_type = ix_type.toIntVect();
443 domain.convert(ix_type);
444 const auto&
dom_hi = ubound(domain);
445 const auto&
dom_lo = lbound(domain);
447 int set_width_x = (iv_type[0]) ? set_width : set_width-1;
448 int set_width_y = (iv_type[1]) ? set_width : set_width-1;
451 #pragma omp parallel if (Gpu::notInLaunchRegion())
453 for (MFIter mfi(S_old_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
454 Box tbx = mfi.tilebox();
455 Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
460 Array4<Real> rhs_arr; Array4<Real> data_arr;
461 Array4<Real> arr_xlo; Array4<Real> arr_xhi;
462 Array4<Real> arr_ylo; Array4<Real> arr_yhi;
464 arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
465 arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
468 }
else if (ivar == ivarV) {
469 arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
470 arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
473 }
else if (ivar == ivarT){
474 arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
475 arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
482 if (ivar == ivarT) {
continue; }
485 width, set_width_x, set_width_y,
487 tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi,
488 arr_xlo, arr_xhi, arr_ylo, arr_yhi,
500 if (width > set_width) {
501 auto dx = geom.CellSizeArray();
502 auto ProbLo = geom.ProbLoArray();
503 auto ProbHi = geom.ProbHiArray();
504 for (
int ivar(ivarU); ivar < BdyEnd; ivar++) {
505 int ivar_idx = ivar_map[ivar];
506 int icomp = comp_map[ivar];
508 Box domain = geom.Domain();
509 domain.convert(S_cur_data[ivar_idx].boxArray().ixType());
510 IntVect ng_vect = S_cur_data[ivar_idx].nGrowVect();
513 #pragma omp parallel if (Gpu::notInLaunchRegion())
515 for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
516 Box tbx = mfi.tilebox();
517 Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
523 Array4<Real> rhs_arr; Array4<Real> data_arr;
524 Array4<Real> arr_xlo; Array4<Real> arr_xhi;
525 Array4<Real> arr_ylo; Array4<Real> arr_yhi;
527 arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
528 arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
531 }
else if (ivar == ivarV) {
532 arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
533 arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
536 }
else if (ivar == ivarT){
537 arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
538 arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
545 if (ivar == ivarT) {
continue; }
548 width, dx, ProbLo, ProbHi, F1, F2,
549 tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi,
550 arr_xlo, arr_xhi, arr_ylo, arr_yhi,
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
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 int &set_width, 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
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_set_rhs_in_spec_region(const amrex::Real &dt, const int &icomp, const int &num_var, const int &width, const int &set_width_x, const int &set_width_y, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, 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 > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr)
Definition: ERF_Utils.H:169
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_laplacian_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::Real &F2, 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 > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr)
Definition: ERF_Utils.H:263
@ 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