ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InteriorGhostCells.cpp File Reference
#include <ERF_Utils.H>
Include dependency graph for ERF_InteriorGhostCells.cpp:

Functions

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)
 
void realbdy_bc_bxs_xy (const Box &bx, const Box &domain, const int &set_width, Box &bx_xlo, Box &bx_xhi, Box &bx_ylo, Box &bx_yhi, const IntVect &ng_vect)
 
void realbdy_compute_interior_ghost_rhs (const Real &bdy_time_interval, const Real &time, const Real &delta_t, const Real &stop_time_elapsed, int width, int set_width, bool do_upwind, const Geometry &geom, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_old_data, Vector< MultiFab > &S_cur_data, Vector< Vector< FArrayBox >> &bdy_data_xlo, Vector< Vector< FArrayBox >> &bdy_data_xhi, Vector< Vector< FArrayBox >> &bdy_data_ylo, Vector< Vector< FArrayBox >> &bdy_data_yhi)
 
void fine_compute_interior_ghost_rhs (const Real &time, const Real &delta_t, const int &width, const int &set_width, const Geometry &geom, ERFFillPatcher *FPr_c, ERFFillPatcher *FPr_u, ERFFillPatcher *FPr_v, ERFFillPatcher *FPr_w, Vector< BCRec > &domain_bcs_type, Vector< MultiFab > &S_rhs_f, Vector< MultiFab > &S_data_f)
 

Variables

PhysBCFunctNoOp void_bc
 

Function Documentation

◆ fine_compute_interior_ghost_rhs()

void fine_compute_interior_ghost_rhs ( const Real time,
const Real delta_t,
const int &  width,
const int &  set_width,
const Geometry &  geom,
ERFFillPatcher FPr_c,
ERFFillPatcher FPr_u,
ERFFillPatcher FPr_v,
ERFFillPatcher FPr_w,
Vector< BCRec > &  domain_bcs_type,
Vector< MultiFab > &  S_rhs_f,
Vector< MultiFab > &  S_data_f 
)

Compute the RHS in the fine relaxation zone

Parameters
[in]timecurrent time
[in]delta_ttimestep
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[in]FPr_ccons fine patch container
[in]FPr_uuvel fine patch container
[in]FPr_vvvel fine patch container
[in]FPr_wwvel fine patch container
[in]boxes_at_levelboxes at current level
[in]domain_bcs_typeboundary condition types
[out]S_rhsRHS to be computed here
[in]S_datacurrent value of the solution
615 {
616  BL_PROFILE_REGION("fine_compute_interior_ghost_RHS()");
617 
618  // Relaxation constants
619  Real F1 = 1./(10.*delta_t);
620  Real F2 = 1./(50.*delta_t);
621 
622  // Vector of MFs to hold data (dm differs w/ fine patch)
623  Vector<MultiFab> fmf_p_v;
624 
625  // Loop over the variables
626  for (int ivar_idx = 0; ivar_idx < IntVars::NumTypes; ++ivar_idx)
627  {
628  // Fine mfs
629  MultiFab& fmf = S_data_f[ivar_idx];
630  MultiFab& rhs = S_rhs_f [ivar_idx];
631 
632  // NOTE: These temporary MFs and copy operations are horrible
633  // for memory usage and efficiency. However, we need to
634  // have access to ghost cells in the cons array to convert
635  // from primitive u/v/w to momentum. Furthermore, the BA
636  // for the fine patches in ERFFillPatcher don't match the
637  // BA for the data/RHS. For this reason, the data is copied
638  // to a vector of MFs (with ghost cells) so the BAs match
639  // the BA of data/RHS and we have access to rho to convert
640  // prim to conserved.
641 
642  // Temp MF on box (distribution map differs w/ fine patch)
643  int num_var = fmf.nComp();
644  fmf_p_v.emplace_back(fmf.boxArray(), fmf.DistributionMap(), num_var, fmf.nGrowVect());
645  MultiFab& fmf_p = fmf_p_v[ivar_idx];
646  MultiFab::Copy(fmf_p,fmf, 0, 0, num_var, fmf.nGrowVect());
647 
648  // Integer mask MF
649  int set_mask_val;
650  int relax_mask_val;
651  iMultiFab* mask;
652 
653  // Fill fine patch on interior halo region
654  //==========================================================
655  if (ivar_idx == IntVars::cons)
656  {
657  FPr_c->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
658  mask = FPr_c->GetMask();
659  set_mask_val = FPr_c->GetSetMaskVal();
660  relax_mask_val = FPr_c->GetRelaxMaskVal();
661  }
662  else if (ivar_idx == IntVars::xmom)
663  {
664  FPr_u->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
665  mask = FPr_u->GetMask();
666  set_mask_val = FPr_u->GetSetMaskVal();
667  relax_mask_val = FPr_u->GetRelaxMaskVal();
668 
669 #ifdef _OPENMP
670 #pragma omp parallel if (Gpu::notInLaunchRegion())
671 #endif
672  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
673  {
674  Box tbx = mfi.tilebox();
675 
676  const Array4<Real>& prim_arr = fmf_p.array(mfi);
677  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
678  const Array4<const int>& mask_arr = mask->const_array(mfi);
679 
680  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
681  {
682  if (mask_arr(i,j,k) == relax_mask_val) {
683  Real rho_interp = 0.5 * ( rho_arr(i-1,j,k) + rho_arr(i,j,k) );
684  prim_arr(i,j,k) *= rho_interp;
685  }
686  });
687  } // mfi
688  }
689  else if (ivar_idx == IntVars::ymom)
690  {
691  FPr_v->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
692  mask = FPr_v->GetMask();
693  set_mask_val = FPr_v->GetSetMaskVal();
694  relax_mask_val = FPr_v->GetRelaxMaskVal();
695 
696 #ifdef _OPENMP
697 #pragma omp parallel if (Gpu::notInLaunchRegion())
698 #endif
699  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
700  {
701  Box tbx = mfi.tilebox();
702 
703  const Array4<Real>& prim_arr = fmf_p.array(mfi);
704  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
705  const Array4<const int>& mask_arr = mask->const_array(mfi);
706 
707  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
708  {
709  if (mask_arr(i,j,k) == relax_mask_val) {
710  Real rho_interp = 0.5 * ( rho_arr(i,j-1,k) + rho_arr(i,j,k) );
711  prim_arr(i,j,k) *= rho_interp;
712  }
713  });
714  } // mfi
715  }
716  else if (ivar_idx == IntVars::zmom)
717  {
718  FPr_w->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
719  mask = FPr_w->GetMask();
720  set_mask_val = FPr_w->GetSetMaskVal();
721  relax_mask_val = FPr_w->GetRelaxMaskVal();
722 
723 #ifdef _OPENMP
724 #pragma omp parallel if (Gpu::notInLaunchRegion())
725 #endif
726  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
727  {
728  Box tbx = mfi.tilebox();
729 
730  const Array4<Real>& prim_arr = fmf_p.array(mfi);
731  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
732  const Array4<const int>& mask_arr = mask->const_array(mfi);
733 
734  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
735  {
736  if (mask_arr(i,j,k) == relax_mask_val) {
737  Real rho_interp = 0.5 * ( rho_arr(i,j,k-1) + rho_arr(i,j,k) );
738  prim_arr(i,j,k) *= rho_interp;
739  }
740  });
741  } // mfi
742  } else {
743  Abort("Dont recognize this variable type in fine_compute_interior_ghost_RHS");
744  }
745 
746 
747  // Zero RHS in set region
748  //==========================================================
749 #ifdef _OPENMP
750 #pragma omp parallel if (Gpu::notInLaunchRegion())
751 #endif
752  for ( MFIter mfi(rhs,TilingIfNotGPU()); mfi.isValid(); ++mfi)
753  {
754  Box tbx = mfi.tilebox();
755  const Array4<Real>& rhs_arr = rhs.array(mfi);
756  const Array4<const int>& mask_arr = mask->const_array(mfi);
757 
758  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
759  {
760  if (mask_arr(i,j,k) == set_mask_val) {
761  rhs_arr(i,j,k) = 0.0;
762  }
763  });
764  } // mfi
765 
766  // For Laplacian stencil
767  rhs.FillBoundary(geom.periodicity());
768 
769 
770  // Compute RHS in relaxation region
771  //==========================================================
772 #ifdef _OPENMP
773 #pragma omp parallel if (Gpu::notInLaunchRegion())
774 #endif
775  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
776  {
777  Box tbx = mfi.tilebox();
778  const Array4<Real>& rhs_arr = rhs.array(mfi);
779  const Array4<const Real>& fine_arr = fmf_p.const_array(mfi);
780  const Array4<const Real>& data_arr = fmf.const_array(mfi);
781  const Array4<const int>& mask_arr = mask->const_array(mfi);
782 
783  Box vbx = mfi.validbox();
784  const auto& vbx_lo = lbound(vbx);
785  const auto& vbx_hi = ubound(vbx);
786 
787  int icomp = 0;
788 
789  int Spec_z = set_width;
790  int Relax_z = width - Spec_z;
791  Real num = Real(Spec_z + Relax_z);
792  Real denom = Real(Relax_z - 1);
793  ParallelFor(tbx, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
794  {
795  if (mask_arr(i,j,k) == relax_mask_val) {
796 
797  // Indices
798  Real n_ind(-1); // Set to -1 to quiet compiler warning
799  int ii(width-1); int jj(width-1);
800  bool near_x_lo_wall(false); bool near_x_hi_wall(false);
801  bool near_y_lo_wall(false); bool near_y_hi_wall(false);
802  bool mask_x_found(false); bool mask_y_found(false);
803 
804  // Near x-wall
805  if ((i-vbx_lo.x) < width) {
806  near_x_lo_wall = true;
807  ii = i-vbx_lo.x;
808  if (mask_arr(vbx_lo.x,j,k) == 2) mask_x_found = true;
809  } else if ((vbx_hi.x-i) < width) {
810  near_x_hi_wall = true;
811  ii = vbx_hi.x-i;
812  if (mask_arr(vbx_hi.x,j,k) == 2) mask_x_found = true;
813  }
814 
815  // Near y-wall
816  if ((j-vbx_lo.y) < width) {
817  near_y_lo_wall = true;
818  jj = j-vbx_lo.y;
819  if (mask_arr(i,vbx_lo.y,k) == 2) mask_y_found = true;
820  } else if ((vbx_hi.y-j) < width) {
821  near_y_hi_wall = true;
822  jj = vbx_hi.y-j;
823  if (mask_arr(i,vbx_hi.y,k) == 2) mask_y_found = true;
824  }
825 
826  // Found a nearby masked cell (valid n_ind)
827  if (mask_x_found && mask_y_found) {
828  n_ind = std::min(ii,jj) + 1.0;
829  } else if (mask_x_found) {
830  n_ind = ii + 1.0;
831  } else if (mask_y_found) {
832  n_ind = jj + 1.0;
833  // Pesky corner cell
834  } else {
835  if (near_x_lo_wall || near_x_hi_wall) {
836  Real dj_min{width-1.0};
837  int j_lb = std::max(vbx_lo.y,j-width);
838  int j_ub = std::min(vbx_hi.y,j+width);
839  int li = (near_x_lo_wall) ? vbx_lo.x : vbx_hi.x;
840  for (int lj(j_lb); lj<=j_ub; ++lj) {
841  if (mask_arr(li,lj,k) == 2) {
842  mask_y_found = true;
843  dj_min = std::min(dj_min,(Real) std::abs(lj-j));
844  }
845  }
846  if (mask_y_found) {
847  Real mag = sqrt( Real(dj_min*dj_min + ii*ii) );
848  n_ind = std::min(mag,width-1.0) + 1.0;
849  } else {
850  Abort("Mask not found near x wall!");
851  }
852  } else if (near_y_lo_wall || near_y_hi_wall) {
853  Real di_min{width-1.0};
854  int i_lb = std::max(vbx_lo.x,i-width);
855  int i_ub = std::min(vbx_hi.x,i+width);
856  int lj = (near_y_lo_wall) ? vbx_lo.y : vbx_hi.y;
857  for (int li(i_lb); li<=i_ub; ++li) {
858  if (mask_arr(li,lj,k) == 2) {
859  mask_x_found = true;
860  di_min = std::min(di_min,(Real) std::abs(li-i));
861  }
862  }
863  if (mask_x_found) {
864  Real mag = sqrt( Real(di_min*di_min + jj*jj) );
865  n_ind = std::min(mag,width-1.0) + 1.0;
866  } else {
867  Abort("Mask not found near y wall!");
868  }
869  } else {
870  Abort("Relaxation cell must be near a wall!");
871  }
872  }
873 
874  Real Factor = (num - n_ind)/denom;
875  Real d = data_arr(i ,j ,k ,n+icomp) + delta_t*rhs_arr(i , j , k ,n+icomp);
876  Real d_ip1 = data_arr(i+1,j ,k ,n+icomp) + delta_t*rhs_arr(i+1, j , k ,n+icomp);
877  Real d_im1 = data_arr(i-1,j ,k ,n+icomp) + delta_t*rhs_arr(i-1, j , k ,n+icomp);
878  Real d_jp1 = data_arr(i ,j+1,k ,n+icomp) + delta_t*rhs_arr(i , j+1, k ,n+icomp);
879  Real d_jm1 = data_arr(i ,j-1,k ,n+icomp) + delta_t*rhs_arr(i , j-1, k ,n+icomp);
880  Real delta = fine_arr(i ,j ,k,n) - d;
881  Real delta_xp = fine_arr(i+1,j ,k,n) - d_ip1;
882  Real delta_xm = fine_arr(i-1,j ,k,n) - d_im1;
883  Real delta_yp = fine_arr(i ,j+1,k,n) - d_jp1;
884  Real delta_ym = fine_arr(i ,j-1,k,n) - d_jm1;
885  Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
886  rhs_arr(i,j,k,n) += (F1*delta - F2*Laplacian) * Factor;
887  }
888  });
889  } // mfi
890  } // ivar_idx
891 }
@ num
Definition: ERF_DataStruct.H:24
PhysBCFunctNoOp void_bc
Definition: ERF_InteriorGhostCells.cpp:5
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::iMultiFab * GetMask()
Definition: ERF_FillPatcher.H:43
void FillRelax(amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
Definition: ERF_FillPatcher.H:105
int GetSetMaskVal()
Definition: ERF_FillPatcher.H:39
int GetRelaxMaskVal()
Definition: ERF_FillPatcher.H:41
@ NumTypes
Definition: ERF_IndexDefines.H:162
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
Here is the call graph for this function:

◆ realbdy_bc_bxs_xy()

void realbdy_bc_bxs_xy ( const Box &  bx,
const Box &  domain,
const int &  set_width,
Box &  bx_xlo,
Box &  bx_xhi,
Box &  bx_ylo,
Box &  bx_yhi,
const IntVect &  ng_vect 
)

Get the boxes for looping over interior/exterior ghost cells for use by fillpatch, erf_slow_rhs_pre, and erf_slow_rhs_post.

Parameters
[in]bxbox to intersect with 4 halo regions
[in]domainbox of the whole domain
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[out]bx_xlohalo box at x_lo boundary
[out]bx_xhihalo box at x_hi boundary
[out]bx_ylohalo box at y_lo boundary
[out]bx_yhihalo box at y_hi boundary
[in]ng_vectnumber of ghost cells in each direction
[in]get_int_ngflag to get ghost cells inside the domain
117 {
118  AMREX_ALWAYS_ASSERT(bx.ixType() == domain.ixType());
119 
120  // Domain bounds without ghost cells
121  const auto& dom_lo = lbound(domain);
122  const auto& dom_hi = ubound(domain);
123 
124  // Four boxes matching the domain
125  Box gdom_xlo(domain); Box gdom_xhi(domain);
126  Box gdom_ylo(domain); Box gdom_yhi(domain);
127 
128  // Get offsets from box index type
129  IntVect iv_type = bx.ixType().toIntVect();
130  int offx = (iv_type[0]==1) ? 0 : -1;
131  int offy = (iv_type[1]==1) ? 0 : -1;
132 
133  // Stagger the boxes based upon index type
134  gdom_xlo += IntVect(offx,0,0); gdom_xhi += IntVect(-offx,0,0);
135  gdom_ylo += IntVect(0,offy,0); gdom_yhi += IntVect(0,-offy,0);
136 
137  // Trim the boxes to only include internal ghost cells
138  gdom_xlo.setBig(0,dom_lo.x+set_width+offx-1); gdom_xhi.setSmall(0,dom_hi.x-set_width-offx+1);
139  gdom_ylo.setBig(1,dom_lo.y+set_width+offy-1); gdom_yhi.setSmall(1,dom_hi.y-set_width-offy+1);
140 
141  // Remove overlapping corners from y-face boxes
142  gdom_ylo.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_ylo.setBig(0,gdom_xhi.smallEnd(0)-1);
143  gdom_yhi.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_yhi.setBig(0,gdom_xhi.smallEnd(0)-1);
144 
145  // Grow boxes to get external ghost cells only
146  gdom_xlo.growLo(0,ng_vect[0]+offx); gdom_xhi.growHi(0,ng_vect[0]+offx);
147  gdom_xlo.grow (1,ng_vect[1] ); gdom_xhi.grow (1,ng_vect[1] );
148  gdom_ylo.growLo(1,ng_vect[1]+offy); gdom_yhi.growHi(1,ng_vect[1]+offy);
149 
150  // Populate everything
151  bx_xlo = (bx & gdom_xlo);
152  bx_xhi = (bx & gdom_xhi);
153  bx_ylo = (bx & gdom_ylo);
154  bx_yhi = (bx & gdom_yhi);
155 }
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1

Referenced by realbdy_interior_bxs_xy().

Here is the caller graph for this function:

◆ realbdy_compute_interior_ghost_rhs()

void realbdy_compute_interior_ghost_rhs ( const Real bdy_time_interval,
const Real time,
const Real delta_t,
const Real stop_time_elapsed,
int  width,
int  set_width,
bool  do_upwind,
const Geometry &  geom,
Vector< MultiFab > &  S_rhs,
Vector< MultiFab > &  S_old_data,
Vector< MultiFab > &  S_cur_data,
Vector< Vector< FArrayBox >> &  bdy_data_xlo,
Vector< Vector< FArrayBox >> &  bdy_data_xhi,
Vector< Vector< FArrayBox >> &  bdy_data_ylo,
Vector< Vector< FArrayBox >> &  bdy_data_yhi 
)

Compute the RHS in the relaxation zone

Parameters
[in]bdy_time_intervaltime interval between boundary condition time stamps
[in]timecurrent time
[in]delta_ttimestep
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[in]geomcontainer for geometric information
[out]S_rhsRHS to be computed here
[in]S_datacurrent value of the solution
[in]bdy_data_xloboundary data on interior of low x-face
[in]bdy_data_xhiboundary data on interior of high x-face
[in]bdy_data_yloboundary data on interior of low y-face
[in]bdy_data_yhiboundary data on interior of high y-face
189 {
190  BL_PROFILE_REGION("realbdy_compute_interior_ghost_RHS()");
191 
192  // Relaxation constants
193  Real F1 = 1./delta_t;
194 
195  // Time interpolation
196  Real dT = bdy_time_interval;
197 
198  // NOTE: this is because we define "time" to be time since start_bdy_time
199  Real time_since_start = time;
200 
201  int n_time = static_cast<int>( time_since_start / dT);
202  Real alpha = (time_since_start - n_time * dT) / dT;
203  AMREX_ALWAYS_ASSERT( alpha >= 0. && alpha <= 1.0);
204  Real oma = 1.0 - alpha;
205 
206  int n_time_p1 = n_time + 1;
207  if ((time == stop_time_elapsed) && (alpha==0)) {
208  // stop time coincides with final bdy snapshot -- don't try to read in
209  // another snapshot
210  n_time_p1 = n_time;
211  }
212 
213  /*
214  // UNIT TEST DEBUG
215  oma = 1.0; alpha = 0.0;
216  */
217 
218  // Temporary FABs for storage (owned/filled on all ranks)
219  FArrayBox U_xlo, U_xhi, U_ylo, U_yhi;
220  FArrayBox V_xlo, V_xhi, V_ylo, V_yhi;
221  FArrayBox T_xlo, T_xhi, T_ylo, T_yhi;
222 
223  // Variable index map (WRFBdyVars -> Vars)
224  Vector<int> var_map = {Vars::xvel, Vars::yvel, Vars::cons, Vars::cons};
225  Vector<int> ivar_map = {IntVars::xmom, IntVars::ymom, IntVars::cons, IntVars::cons};
226 
227  // Variable icomp map
228  Vector<int> comp_map = {0, 0, RhoTheta_comp};
229 
230  // Indices
231  int ivarU = RealBdyVars::U;
232  int ivarV = RealBdyVars::V;
233  int ivarT = RealBdyVars::T;
234  int BdyEnd = RealBdyVars::NumTypes-1;
235 
236 
237  // NOTE: The sizing of the temporary BDY FABS is
238  // GLOBAL and occurs over the entire BDY region.
239 
240  // Size the FABs
241  //==========================================================
242  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
243  int ivar_idx = var_map[ivar];
244  Box domain = geom.Domain();
245  auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
246  domain.convert(ixtype);
247 
248  // NOTE: Ghost cells needed for idx type mismatch between mask and data (do_upwind)
249  IntVect ng_vect(1,1,0);
250  Box gdom(domain); gdom.grow(ng_vect);
251  Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
252  realbdy_interior_bxs_xy(gdom, domain, width,
253  bx_xlo, bx_xhi,
254  bx_ylo, bx_yhi,
255  0, ng_vect, true);
256 
257  // Size the FABs
258  if (ivar == ivarU) {
259  U_xlo.resize(bx_xlo,1,The_Async_Arena()); U_xhi.resize(bx_xhi,1,The_Async_Arena());
260  U_ylo.resize(bx_ylo,1,The_Async_Arena()); U_yhi.resize(bx_yhi,1,The_Async_Arena());
261  } else if (ivar == ivarV) {
262  V_xlo.resize(bx_xlo,1,The_Async_Arena()); V_xhi.resize(bx_xhi,1,The_Async_Arena());
263  V_ylo.resize(bx_ylo,1,The_Async_Arena()); V_yhi.resize(bx_yhi,1,The_Async_Arena());
264  } else if (ivar == ivarT){
265  T_xlo.resize(bx_xlo,1,The_Async_Arena()); T_xhi.resize(bx_xhi,1,The_Async_Arena());
266  T_ylo.resize(bx_ylo,1,The_Async_Arena()); T_yhi.resize(bx_yhi,1,The_Async_Arena());
267  } else {
268  continue;
269  }
270  } // ivar
271 
272 
273  // NOTE: These operations use the BDY FABS and RHO. The
274  // use of RHO to go from PRIM -> CONS requires that
275  // these operations be LOCAL. So we have allocated
276  // enough space to do global operations (1 rank) but
277  // will fill a subset of that data that the rank owns.
278 
279  // Populate FABs from bdy interpolation (primitive vars)
280  //==========================================================
281  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
282  int ivar_idx = var_map[ivar];
283  Box domain = geom.Domain();
284  auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
285  domain.convert(ixtype);
286  const auto& dom_lo = lbound(domain);
287  const auto& dom_hi = ubound(domain);
288 
289 #ifdef _OPENMP
290 #pragma omp parallel if (Gpu::notInLaunchRegion())
291 #endif
292  for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
293  // NOTE: Ghost cells needed for idx type mismatch between mask and data (do_upwind)
294  IntVect ng_vect(1,1,0);
295  Box gtbx = grow(mfi.tilebox(ixtype.toIntVect()),ng_vect);
296  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
297  realbdy_interior_bxs_xy(gtbx, domain, width,
298  tbx_xlo, tbx_xhi,
299  tbx_ylo, tbx_yhi,
300  0, ng_vect, true);
301 
302  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
303  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
304  if (ivar == ivarU) {
305  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
306  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
307  } else if (ivar == ivarV) {
308  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
309  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
310  } else if (ivar == ivarT){
311  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
312  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
313  } else {
314  continue;
315  }
316 
317  // Boundary data at fixed time intervals
318  const auto& bdatxlo_n = bdy_data_xlo[n_time ][ivar].const_array();
319  const auto& bdatxlo_np1 = bdy_data_xlo[n_time_p1][ivar].const_array();
320  const auto& bdatxhi_n = bdy_data_xhi[n_time ][ivar].const_array();
321  const auto& bdatxhi_np1 = bdy_data_xhi[n_time_p1][ivar].const_array();
322  const auto& bdatylo_n = bdy_data_ylo[n_time ][ivar].const_array();
323  const auto& bdatylo_np1 = bdy_data_ylo[n_time_p1][ivar].const_array();
324  const auto& bdatyhi_n = bdy_data_yhi[n_time ][ivar].const_array();
325  const auto& bdatyhi_np1 = bdy_data_yhi[n_time_p1][ivar].const_array();
326 
327  // Current density to convert to conserved vars
328  Array4<Real> r_arr = S_cur_data[IntVars::cons].array(mfi);
329 
330  // Limiting offset
331  int offset = set_width - 1;
332  if (width > set_width) offset = width - 1;
333 
334  // Populate with interpolation (protect from ghost cells)
335  ParallelFor(tbx_xlo, tbx_xhi,
336  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
337  {
338  int ii = std::max(i , dom_lo.x);
339  ii = std::min(ii, dom_lo.x+offset);
340  int jj = std::max(j , dom_lo.y);
341  jj = std::min(jj, dom_hi.y);
342 
343  Real rho_interp;
344  if (ivar==ivarU) {
345  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
346  } else if (ivar==ivarV) {
347  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
348  } else {
349  rho_interp = r_arr(i,j,k);
350  }
351  arr_xlo(i,j,k) = rho_interp * ( oma * bdatxlo_n (ii,jj,k,0)
352  + alpha * bdatxlo_np1(ii,jj,k,0) );
353  },
354  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
355  {
356  int ii = std::max(i , dom_hi.x-offset);
357  ii = std::min(ii, dom_hi.x);
358  int jj = std::max(j , dom_lo.y);
359  jj = std::min(jj, dom_hi.y);
360 
361  Real rho_interp;
362  if (ivar==ivarU) {
363  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
364  } else if (ivar==ivarV) {
365  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
366  } else {
367  rho_interp = r_arr(i,j,k);
368  }
369  arr_xhi(i,j,k) = rho_interp * ( oma * bdatxhi_n (ii,jj,k,0)
370  + alpha * bdatxhi_np1(ii,jj,k,0) );
371  });
372 
373  ParallelFor(tbx_ylo, tbx_yhi,
374  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
375  {
376  int ii = std::max(i , dom_lo.x);
377  ii = std::min(ii, dom_hi.x);
378  int jj = std::max(j , dom_lo.y);
379  jj = std::min(jj, dom_lo.y+offset);
380 
381  Real rho_interp;
382  if (ivar==ivarU) {
383  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
384  } else if (ivar==ivarV) {
385  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
386  } else {
387  rho_interp = r_arr(i,j,k);
388  }
389  arr_ylo(i,j,k) = rho_interp * ( oma * bdatylo_n (ii,jj,k,0)
390  + alpha * bdatylo_np1(ii,jj,k,0) );
391  },
392  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
393  {
394  int ii = std::max(i , dom_lo.x);
395  ii = std::min(ii, dom_hi.x);
396  int jj = std::max(j , dom_hi.y-offset);
397  jj = std::min(jj, dom_hi.y);
398 
399  Real rho_interp;
400  if (ivar==ivarU) {
401  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
402  } else if (ivar==ivarV) {
403  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
404  } else {
405  rho_interp = r_arr(i,j,k);
406  }
407  arr_yhi(i,j,k) = rho_interp * ( oma * bdatyhi_n (ii,jj,k,0)
408  + alpha * bdatyhi_np1(ii,jj,k,0) );
409  });
410  } // mfi
411  } // ivar
412 
413 
414  // NOTE: These operations use current RHS, so they are
415  // LOCAL and occur over the data owned by a given rank.
416 
417  // Compute RHS in specified region
418  //==========================================================
419  if (set_width > 0) {
420  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
421  int ivar_idx = ivar_map[ivar];
422  int icomp = comp_map[ivar];
423 
424  Box domain = geom.Domain();
425  auto ix_type = S_old_data[ivar_idx].boxArray().ixType();
426  auto iv_type = ix_type.toIntVect();
427  domain.convert(ix_type);
428 
429  int set_width_x = (iv_type[0]) ? set_width : set_width-1;
430  int set_width_y = (iv_type[1]) ? set_width : set_width-1;
431 
432 #ifdef _OPENMP
433 #pragma omp parallel if (Gpu::notInLaunchRegion())
434 #endif
435  for (MFIter mfi(S_old_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
436  Box tbx = mfi.tilebox();
437  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
438  realbdy_interior_bxs_xy(tbx, domain, width,
439  tbx_xlo, tbx_xhi,
440  tbx_ylo, tbx_yhi);
441 
442  Array4<Real> rhs_arr; Array4<Real> data_arr;
443  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
444  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
445  if (ivar == ivarU) {
446  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
447  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
448  rhs_arr = S_rhs[IntVars::xmom].array(mfi);
449  data_arr = S_old_data[IntVars::xmom].array(mfi);
450  } else if (ivar == ivarV) {
451  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
452  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
453  rhs_arr = S_rhs[IntVars::ymom].array(mfi);
454  data_arr = S_old_data[IntVars::ymom].array(mfi);
455  } else if (ivar == ivarT){
456  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
457  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
458  rhs_arr = S_rhs[IntVars::cons].array(mfi);
459  data_arr = S_old_data[IntVars::cons].array(mfi);
460  } else {
461  continue;
462  }
463 
464  Array4<Real> u_xlo = U_xlo.array(); Array4<Real> u_xhi = U_xhi.array();
465  Array4<Real> v_xlo = V_xlo.array(); Array4<Real> v_xhi = V_xhi.array();
466  Array4<Real> v_ylo = V_ylo.array(); Array4<Real> v_yhi = V_yhi.array();
467 
468  realbdy_set_rhs_in_spec_region(delta_t, icomp, 1,
469  width, set_width_x, set_width_y,
470  domain, geom.Domain(),
471  tbx_xlo , tbx_xhi , tbx_ylo , tbx_yhi ,
472  arr_xlo , arr_xhi , arr_ylo , arr_yhi ,
473  u_xlo, u_xhi, v_xlo, v_xhi, v_ylo, v_yhi,
474  data_arr, rhs_arr, do_upwind);
475 
476  } // mfi
477  } // ivar
478  } // set_width
479 
480  // NOTE: These operations use current density, so they are
481  // LOCAL and occur over the data owned by a given rank
482 
483  // Compute RHS in relaxation region
484  //==========================================================
485  if (width > set_width) {
486  auto dx = geom.CellSizeArray();
487  auto ProbLo = geom.ProbLoArray();
488  auto ProbHi = geom.ProbHiArray();
489  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
490  int ivar_idx = ivar_map[ivar];
491  int icomp = comp_map[ivar];
492 
493  Box domain = geom.Domain();
494  domain.convert(S_cur_data[ivar_idx].boxArray().ixType());
495  IntVect ng_vect(0);
496 
497 #ifdef _OPENMP
498 #pragma omp parallel if (Gpu::notInLaunchRegion())
499 #endif
500  for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
501  Box tbx = mfi.tilebox();
502  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
503  realbdy_interior_bxs_xy(tbx, domain, width,
504  tbx_xlo, tbx_xhi,
505  tbx_ylo, tbx_yhi,
506  set_width, ng_vect);
507 
508  Array4<Real> rhs_arr; Array4<Real> data_arr;
509  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
510  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
511  if (ivar == ivarU) {
512  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
513  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
514  rhs_arr = S_rhs[IntVars::xmom].array(mfi);
515  data_arr = S_cur_data[IntVars::xmom].array(mfi);
516  } else if (ivar == ivarV) {
517  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
518  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
519  rhs_arr = S_rhs[IntVars::ymom].array(mfi);
520  data_arr = S_cur_data[IntVars::ymom].array(mfi);
521  } else if (ivar == ivarT){
522  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
523  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
524  rhs_arr = S_rhs[IntVars::cons].array(mfi);
525  data_arr = S_cur_data[IntVars::cons].array(mfi);
526  } else {
527  continue;
528  }
529 
530  Array4<Real> u_xlo = U_xlo.array(); Array4<Real> u_xhi = U_xhi.array();
531  Array4<Real> v_xlo = V_xlo.array(); Array4<Real> v_xhi = V_xhi.array();
532  Array4<Real> v_ylo = V_ylo.array(); Array4<Real> v_yhi = V_yhi.array();
533 
535  width, dx, ProbLo, ProbHi, F1, geom.Domain(),
536  tbx_xlo , tbx_xhi , tbx_ylo , tbx_yhi ,
537  arr_xlo , arr_xhi , arr_ylo , arr_yhi ,
538  u_xlo, u_xhi, v_xlo, v_xhi, v_ylo, v_yhi,
539  data_arr, rhs_arr, do_upwind);
540 
541  /*
542  // UNIT TEST DEBUG
543  realbdy_interior_bxs_xy(tbx, domain, width,
544  tbx_xlo, tbx_xhi,
545  tbx_ylo, tbx_yhi);
546  ParallelFor(tbx_xlo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
547  {
548  if (std::fabs(arr_xlo(i,j,k) - data_arr(i,j,k,icomp)) > 0.01) {
549  Print() << "ERROR XLO: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
550  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_xlo(i,j,k) << "\n";
551  exit(0);
552  }
553  });
554  ParallelFor(tbx_xhi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
555  {
556  if (std::fabs(arr_xhi(i,j,k) - data_arr(i,j,k,icomp)) > 0.01) {
557  Print() << "ERROR XHI: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
558  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_xhi(i,j,k) << "\n";
559  exit(0);
560  }
561  });
562  ParallelFor(tbx_ylo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
563  {
564  if (std::fabs(arr_ylo(i,j,k) - data_arr(i,j,k,icomp)) > 0.01) {
565  Print() << "ERROR YLO: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
566  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_ylo(i,j,k) << "\n";
567  exit(0);
568  }
569  });
570  ParallelFor(tbx_yhi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
571  {
572  if (std::fabs(arr_yhi(i,j,k)-data_arr(i,j,k,icomp)) > 0.01) {
573  Print() << "ERROR YHI: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
574  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_yhi(i,j,k) << "\n";
575  exit(0);
576  }
577  });
578  */
579  } // mfi
580  } // ivar
581  } // width
582  //ParallelDescriptor::Barrier();
583  //exit(0);
584 }
#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_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:323
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::Box &domain, 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:183
@ 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
Here is the call graph for this function:

◆ realbdy_interior_bxs_xy()

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 
)

Get the boxes for looping over interior/exterior ghost cells for use by fillpatch, erf_slow_rhs_pre, and erf_slow_rhs_post.

Parameters
[in]bxbox to intersect with 4 halo regions
[in]domainbox of the whole domain
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[out]bx_xlohalo box at x_lo boundary
[out]bx_xhihalo box at x_hi boundary
[out]bx_ylohalo box at y_lo boundary
[out]bx_yhihalo box at y_hi boundary
[in]ng_vectnumber of ghost cells in each direction
[in]get_int_ngflag to get ghost cells inside the domain
33 {
34  AMREX_ALWAYS_ASSERT(bx.ixType() == domain.ixType());
35 
36  //==================================================================
37  // NOTE: X-face boxes take ownership of the overlapping region.
38  // With exterior ghost cells (ng_vect != 0), the x-face
39  // boxes will have exterior ghost cells in both x & y.
40  //==================================================================
41 
42  // Domain bounds without ghost cells
43  const auto& dom_lo = lbound(domain);
44  const auto& dom_hi = ubound(domain);
45 
46  // Four boxes matching the domain
47  Box gdom_xlo(domain); Box gdom_xhi(domain);
48  Box gdom_ylo(domain); Box gdom_yhi(domain);
49 
50  // Trim the boxes to only include internal ghost cells
51  gdom_xlo.setBig(0,dom_lo.x+width-1); gdom_xhi.setSmall(0,dom_hi.x-width+1);
52  gdom_ylo.setBig(1,dom_lo.y+width-1); gdom_yhi.setSmall(1,dom_hi.y-width+1);
53 
54  // Remove overlapping corners from y-face boxes
55  gdom_ylo.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_ylo.setBig(0,gdom_xhi.smallEnd(0)-1);
56  gdom_yhi.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_yhi.setBig(0,gdom_xhi.smallEnd(0)-1);
57 
58  if (set_width>0) {
59  // Cut out the set region (uses offsets)
60  IntVect iv_type = bx.ixType().toIntVect();
61  Box sdom_xlo,sdom_xhi,sdom_ylo,sdom_yhi;
62  realbdy_bc_bxs_xy (grow(domain,ng_vect), domain, set_width,
63  sdom_xlo, sdom_xhi,
64  sdom_ylo, sdom_yhi,
65  ng_vect);
66  gdom_xlo.setSmall(0,sdom_xlo.bigEnd(0)+1); gdom_xhi.setBig(0,sdom_xhi.smallEnd(0)-1);
67  gdom_ylo.setSmall(1,sdom_ylo.bigEnd(1)+1); gdom_yhi.setBig(1,sdom_yhi.smallEnd(1)-1);
68  if (iv_type[1]==1) {
69  gdom_xlo.setSmall(1,sdom_ylo.bigEnd(1)+1); gdom_xhi.setSmall(1,sdom_ylo.bigEnd(1)+1);
70  gdom_xlo.setBig(1,sdom_yhi.smallEnd(1)-1); gdom_xhi.setBig(1,sdom_yhi.smallEnd(1)-1);
71  }
72  } else {
73  // Grow boxes to get external ghost cells only
74  gdom_xlo.growLo(0,ng_vect[0]); gdom_xhi.growHi(0,ng_vect[0]);
75  gdom_xlo.grow (1,ng_vect[1]); gdom_xhi.grow (1,ng_vect[1]);
76  gdom_ylo.growLo(1,ng_vect[1]); gdom_yhi.growHi(1,ng_vect[1]);
77  }
78 
79  // Grow boxes to get internal ghost cells
80  if (get_int_ng) {
81  gdom_xlo.growHi(0,ng_vect[0]); gdom_xhi.growLo(0,ng_vect[0]);
82  gdom_ylo.grow (0,ng_vect[0]); gdom_yhi.grow (0,ng_vect[0]);
83  gdom_ylo.growHi(1,ng_vect[1]); gdom_yhi.growLo(1,ng_vect[1]);
84  }
85 
86  // Populate everything
87  bx_xlo = (bx & gdom_xlo);
88  bx_xhi = (bx & gdom_xhi);
89  bx_ylo = (bx & gdom_ylo);
90  bx_yhi = (bx & gdom_yhi);
91 }
void realbdy_bc_bxs_xy(const Box &bx, const Box &domain, const int &set_width, Box &bx_xlo, Box &bx_xhi, Box &bx_ylo, Box &bx_yhi, const IntVect &ng_vect)
Definition: ERF_InteriorGhostCells.cpp:109

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ void_bc

PhysBCFunctNoOp void_bc