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 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 &time, const Real &delta_t, const Real &start_bdy_time, const Real &final_bdy_time, const Real &bdy_time_interval, const Real &nudge_factor, int width, const Geometry &geom, Vector< MultiFab > &S_rhs, 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, std::unique_ptr< ReadBndryPlanes > &m_r2d)
 
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 (elapsed) 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
642 {
643  BL_PROFILE_REGION("fine_compute_interior_ghost_RHS()");
644 
645  // Relaxation constants
646  Real F1 = one/(Real(10.)*delta_t);
647  Real F2 = one/(Real(50.)*delta_t);
648 
649  // Vector of MFs to hold data (dm differs w/ fine patch)
650  Vector<MultiFab> fmf_p_v;
651 
652  // Loop over the variables
653  for (int ivar_idx = 0; ivar_idx < IntVars::NumTypes; ++ivar_idx)
654  {
655  // Fine mfs
656  MultiFab& fmf = S_data_f[ivar_idx];
657  MultiFab& rhs = S_rhs_f [ivar_idx];
658 
659  // NOTE: These temporary MFs and copy operations are horrible
660  // for memory usage and efficiency. However, we need to
661  // have access to ghost cells in the cons array to convert
662  // from primitive u/v/w to momentum. Furthermore, the BA
663  // for the fine patches in ERFFillPatcher don't match the
664  // BA for the data/RHS. For this reason, the data is copied
665  // to a vector of MFs (with ghost cells) so the BAs match
666  // the BA of data/RHS and we have access to rho to convert
667  // prim to conserved.
668 
669  // Temp MF on box (distribution map differs w/ fine patch)
670  int num_var = fmf.nComp();
671  fmf_p_v.emplace_back(fmf.boxArray(), fmf.DistributionMap(), num_var, fmf.nGrowVect());
672  MultiFab& fmf_p = fmf_p_v[ivar_idx];
673  MultiFab::Copy(fmf_p,fmf, 0, 0, num_var, fmf.nGrowVect());
674 
675  // Integer mask MF
676  int set_mask_val;
677  int relax_mask_val;
678  iMultiFab* mask;
679 
680  // Fill fine patch on interior halo region
681  //==========================================================
682  if (ivar_idx == IntVars::cons)
683  {
684  FPr_c->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
685  mask = FPr_c->GetMask();
686  set_mask_val = FPr_c->GetSetMaskVal();
687  relax_mask_val = FPr_c->GetRelaxMaskVal();
688  }
689  else if (ivar_idx == IntVars::xmom)
690  {
691  FPr_u->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
692  mask = FPr_u->GetMask();
693  set_mask_val = FPr_u->GetSetMaskVal();
694  relax_mask_val = FPr_u->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 = myhalf * ( rho_arr(i-1,j,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::ymom)
717  {
718  FPr_v->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
719  mask = FPr_v->GetMask();
720  set_mask_val = FPr_v->GetSetMaskVal();
721  relax_mask_val = FPr_v->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 = myhalf * ( rho_arr(i,j-1,k) + rho_arr(i,j,k) );
738  prim_arr(i,j,k) *= rho_interp;
739  }
740  });
741  } // mfi
742  }
743  else if (ivar_idx == IntVars::zmom)
744  {
745  FPr_w->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
746  mask = FPr_w->GetMask();
747  set_mask_val = FPr_w->GetSetMaskVal();
748  relax_mask_val = FPr_w->GetRelaxMaskVal();
749 
750 #ifdef _OPENMP
751 #pragma omp parallel if (Gpu::notInLaunchRegion())
752 #endif
753  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
754  {
755  Box tbx = mfi.tilebox();
756 
757  const Array4<Real>& prim_arr = fmf_p.array(mfi);
758  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
759  const Array4<const int>& mask_arr = mask->const_array(mfi);
760 
761  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
762  {
763  if (mask_arr(i,j,k) == relax_mask_val) {
764  Real rho_interp = myhalf * ( rho_arr(i,j,k-1) + rho_arr(i,j,k) );
765  prim_arr(i,j,k) *= rho_interp;
766  }
767  });
768  } // mfi
769  } else {
770  Abort("Dont recognize this variable type in fine_compute_interior_ghost_RHS");
771  }
772 
773 
774  // Zero RHS in set region
775  //==========================================================
776 #ifdef _OPENMP
777 #pragma omp parallel if (Gpu::notInLaunchRegion())
778 #endif
779  for ( MFIter mfi(rhs,TilingIfNotGPU()); mfi.isValid(); ++mfi)
780  {
781  Box tbx = mfi.tilebox();
782  const Array4<Real>& rhs_arr = rhs.array(mfi);
783  const Array4<const int>& mask_arr = mask->const_array(mfi);
784 
785  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
786  {
787  if (mask_arr(i,j,k) == set_mask_val) {
788  rhs_arr(i,j,k) = zero;
789  }
790  });
791  } // mfi
792 
793  // For Laplacian stencil
794  rhs.FillBoundary(geom.periodicity());
795 
796 
797  // Compute RHS in relaxation region
798  //==========================================================
799 #ifdef _OPENMP
800 #pragma omp parallel if (Gpu::notInLaunchRegion())
801 #endif
802  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
803  {
804  Box tbx = mfi.tilebox();
805  const Array4<Real>& rhs_arr = rhs.array(mfi);
806  const Array4<const Real>& fine_arr = fmf_p.const_array(mfi);
807  const Array4<const Real>& data_arr = fmf.const_array(mfi);
808  const Array4<const int>& mask_arr = mask->const_array(mfi);
809 
810  Box vbx = mfi.validbox();
811  const auto& vbx_lo = lbound(vbx);
812  const auto& vbx_hi = ubound(vbx);
813 
814  int icomp = 0;
815 
816  int Spec_z = set_width;
817  int Relax_z = width - Spec_z;
818  Real num = Real(Spec_z + Relax_z);
819  Real denom = Real(Relax_z - 1);
820  ParallelFor(tbx, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
821  {
822  if (mask_arr(i,j,k) == relax_mask_val) {
823 
824  // Indices
825  Real n_ind(-1); // Set to -1 to quiet compiler warning
826  int ii(width-1); int jj(width-1);
827  bool near_x_lo_wall(false); bool near_x_hi_wall(false);
828  bool near_y_lo_wall(false); bool near_y_hi_wall(false);
829  bool mask_x_found(false); bool mask_y_found(false);
830 
831  // Near x-wall
832  if ((i-vbx_lo.x) < width) {
833  near_x_lo_wall = true;
834  ii = i-vbx_lo.x;
835  if (mask_arr(vbx_lo.x,j,k) == 2) mask_x_found = true;
836  } else if ((vbx_hi.x-i) < width) {
837  near_x_hi_wall = true;
838  ii = vbx_hi.x-i;
839  if (mask_arr(vbx_hi.x,j,k) == 2) mask_x_found = true;
840  }
841 
842  // Near y-wall
843  if ((j-vbx_lo.y) < width) {
844  near_y_lo_wall = true;
845  jj = j-vbx_lo.y;
846  if (mask_arr(i,vbx_lo.y,k) == 2) mask_y_found = true;
847  } else if ((vbx_hi.y-j) < width) {
848  near_y_hi_wall = true;
849  jj = vbx_hi.y-j;
850  if (mask_arr(i,vbx_hi.y,k) == 2) mask_y_found = true;
851  }
852 
853  // Found a nearby masked cell (valid n_ind)
854  if (mask_x_found && mask_y_found) {
855  n_ind = std::min(ii,jj) + one;
856  } else if (mask_x_found) {
857  n_ind = ii + one;
858  } else if (mask_y_found) {
859  n_ind = jj + one;
860  // Pesky corner cell
861  } else {
862  if (near_x_lo_wall || near_x_hi_wall) {
863  Real dj_min{width-one};
864  int j_lb = std::max(vbx_lo.y,j-width);
865  int j_ub = std::min(vbx_hi.y,j+width);
866  int li = (near_x_lo_wall) ? vbx_lo.x : vbx_hi.x;
867  for (int lj(j_lb); lj<=j_ub; ++lj) {
868  if (mask_arr(li,lj,k) == 2) {
869  mask_y_found = true;
870  dj_min = std::min(dj_min,(Real) std::abs(lj-j));
871  }
872  }
873  if (mask_y_found) {
874  Real mag = std::sqrt( Real(dj_min*dj_min + ii*ii) );
875  n_ind = std::min(mag,width-one) + one;
876  } else {
877  Abort("Mask not found near x wall!");
878  }
879  } else if (near_y_lo_wall || near_y_hi_wall) {
880  Real di_min{width-one};
881  int i_lb = std::max(vbx_lo.x,i-width);
882  int i_ub = std::min(vbx_hi.x,i+width);
883  int lj = (near_y_lo_wall) ? vbx_lo.y : vbx_hi.y;
884  for (int li(i_lb); li<=i_ub; ++li) {
885  if (mask_arr(li,lj,k) == 2) {
886  mask_x_found = true;
887  di_min = std::min(di_min,(Real) std::abs(li-i));
888  }
889  }
890  if (mask_x_found) {
891  Real mag = std::sqrt( Real(di_min*di_min + jj*jj) );
892  n_ind = std::min(mag,width-one) + one;
893  } else {
894  Abort("Mask not found near y wall!");
895  }
896  } else {
897  Abort("Relaxation cell must be near a wall!");
898  }
899  }
900 
901  Real Factor = (num - n_ind)/denom;
902  Real d = data_arr(i ,j ,k ,n+icomp) + delta_t*rhs_arr(i , j , k ,n+icomp);
903  Real d_ip1 = data_arr(i+1,j ,k ,n+icomp) + delta_t*rhs_arr(i+1, j , k ,n+icomp);
904  Real d_im1 = data_arr(i-1,j ,k ,n+icomp) + delta_t*rhs_arr(i-1, j , k ,n+icomp);
905  Real d_jp1 = data_arr(i ,j+1,k ,n+icomp) + delta_t*rhs_arr(i , j+1, k ,n+icomp);
906  Real d_jm1 = data_arr(i ,j-1,k ,n+icomp) + delta_t*rhs_arr(i , j-1, k ,n+icomp);
907  Real delta = fine_arr(i ,j ,k,n) - d;
908  Real delta_xp = fine_arr(i+1,j ,k,n) - d_ip1;
909  Real delta_xm = fine_arr(i-1,j ,k,n) - d_im1;
910  Real delta_yp = fine_arr(i ,j+1,k,n) - d_jp1;
911  Real delta_ym = fine_arr(i ,j-1,k,n) - d_jm1;
912  Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - Real(4.0)*delta;
913  rhs_arr(i,j,k,n) += (F1*delta - F2*Laplacian) * Factor;
914  }
915  });
916  } // mfi
917  } // ivar_idx
918 }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
@ num
Definition: ERF_DataStruct.H:24
for(int i=0;i< m_num_species;i++)
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:48
PhysBCFunctNoOp void_bc
Definition: ERF_InteriorGhostCells.cpp:5
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
auto rho_arr
Definition: ERF_UpdateWSubsidence_SineMassFlux.H:3
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:196
@ ymom
Definition: ERF_IndexDefines.H:194
@ cons
Definition: ERF_IndexDefines.H:192
@ zmom
Definition: ERF_IndexDefines.H:195
@ xmom
Definition: ERF_IndexDefines.H:193
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
101 {
102  AMREX_ALWAYS_ASSERT(bx.ixType() == domain.ixType());
103 
104  // Domain bounds without ghost cells
105  const auto& dom_lo = lbound(domain);
106  const auto& dom_hi = ubound(domain);
107 
108  // Four boxes matching the domain
109  Box gdom_xlo(domain); Box gdom_xhi(domain);
110  Box gdom_ylo(domain); Box gdom_yhi(domain);
111 
112  // Get offsets from box index type
113  IntVect iv_type = bx.ixType().toIntVect();
114  int offx = (iv_type[0]==1) ? 0 : -1;
115  int offy = (iv_type[1]==1) ? 0 : -1;
116 
117  // Stagger the boxes based upon index type
118  gdom_xlo += IntVect(offx,0,0); gdom_xhi += IntVect(-offx,0,0);
119  gdom_ylo += IntVect(0,offy,0); gdom_yhi += IntVect(0,-offy,0);
120 
121  // Trim the boxes to only include internal ghost cells
122  gdom_xlo.setBig(0,dom_lo.x+set_width+offx-1); gdom_xhi.setSmall(0,dom_hi.x-set_width-offx+1);
123  gdom_ylo.setBig(1,dom_lo.y+set_width+offy-1); gdom_yhi.setSmall(1,dom_hi.y-set_width-offy+1);
124 
125  // Remove overlapping corners from y-face boxes
126  gdom_ylo.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_ylo.setBig(0,gdom_xhi.smallEnd(0)-1);
127  gdom_yhi.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_yhi.setBig(0,gdom_xhi.smallEnd(0)-1);
128 
129  // Grow boxes to get external ghost cells only
130  gdom_xlo.growLo(0,ng_vect[0]+offx); gdom_xhi.growHi(0,ng_vect[0]+offx);
131  gdom_xlo.grow (1,ng_vect[1] ); gdom_xhi.grow (1,ng_vect[1] );
132  gdom_ylo.growLo(1,ng_vect[1]+offy); gdom_yhi.growHi(1,ng_vect[1]+offy);
133 
134  // Populate everything
135  bx_xlo = (bx & gdom_xlo);
136  bx_xhi = (bx & gdom_xhi);
137  bx_ylo = (bx & gdom_ylo);
138  bx_yhi = (bx & gdom_yhi);
139 }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
Here is the call graph for this function:

◆ realbdy_compute_interior_ghost_rhs()

void realbdy_compute_interior_ghost_rhs ( const Real time,
const Real delta_t,
const Real start_bdy_time,
const Real final_bdy_time,
const Real bdy_time_interval,
const Real nudge_factor,
int  width,
const Geometry &  geom,
Vector< MultiFab > &  S_rhs,
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,
std::unique_ptr< ReadBndryPlanes > &  m_r2d 
)

Compute the RHS in the relaxation zone

Parameters
[in]timecurrent (total) time
[in]delta_ttimestep
[in]start_bdy_timefull time of the first time slice of boundary data
[in]final_bdy_timefull time of the last time slice of boundary data
[in]bdy_time_intervaltime interval between boundary condition time stamps
[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
176 {
177  BL_PROFILE_REGION("realbdy_compute_interior_ghost_RHS()");
178 
179  //
180  // Note that time (= start_time+old_stage_time) is measured as total time
181  // start_bdy_time and final_bdy_time are also measured as total time
182  //
183 
184  // Get bndry data if we have it
185  Vector<int> bnd_map = {BCVars::xvel_bc, BCVars::yvel_bc, BCVars::RhoTheta_bc_comp};
186  Array4<Real> bdatxlo, bdatxhi, bdatylo, bdatyhi;
187  Array4<Real> btenxlo, btenxhi, btenylo, btenyhi;
188  if (m_r2d) {
189  // Index is [plane orientation] and [level]
190  Vector<std::unique_ptr<PlaneVector>>& bndry_data = m_r2d->interp_in_time(time);
191  bdatxlo = (*bndry_data[0])[0].array();
192  bdatylo = (*bndry_data[1])[0].array();
193  bdatxhi = (*bndry_data[3])[0].array();
194  bdatyhi = (*bndry_data[4])[0].array();
195 
196  Vector<std::unique_ptr<PlaneVector>>& bndry_tend = m_r2d->get_tendency(time);
197  btenxlo = (*bndry_tend[0])[0].array();
198  btenylo = (*bndry_tend[1])[0].array();
199  btenxhi = (*bndry_tend[3])[0].array();
200  btenyhi = (*bndry_tend[4])[0].array();
201  }
202 
203  // Relaxation constants
204  Real F1 = one/(nudge_factor*delta_t);
205 
206  // Time interpolation
207  Real dT = bdy_time_interval;
208 
209  int n_time = static_cast<int>( (time-start_bdy_time) / dT);
210  int n_time_p1 = n_time + 1;
211  Real alpha = ((time-start_bdy_time) - n_time * dT) / dT;
212 
213  // Do not over run the last bdy file
214  if (time >= final_bdy_time) {
215  n_time = static_cast<int>( (final_bdy_time - start_bdy_time)/ dT);
216  n_time_p1 = n_time;
217  alpha = zero;
218  }
219 
221  Real oma = one - alpha;
222 
223  // Temporary FABs for storage (owned/filled on all ranks)
224  FArrayBox U_xlo, U_xhi, U_ylo, U_yhi;
225  FArrayBox V_xlo, V_xhi, V_ylo, V_yhi;
226  FArrayBox T_xlo, T_xhi, T_ylo, T_yhi;
227 
228  // Variable index map (WRFBdyVars -> Vars)
229  Vector<int> var_map = {Vars::xvel, Vars::yvel, Vars::cons };
230  Vector<int> ivar_map = {IntVars::xmom, IntVars::ymom, IntVars::cons};
231 
232  // Variable icomp map
233  Vector<int> comp_map = {0, 0, RhoTheta_comp};
234 
235  // Indices
236  int ivarU = RealBdyVars::U;
237  int ivarV = RealBdyVars::V;
238  int ivarT = RealBdyVars::T;
239  int BdyEnd = RealBdyVars::NumTypes-1;
240 
241 
242  // NOTE: The sizing of the temporary BDY FABS is
243  // GLOBAL and occurs over the entire BDY region.
244 
245  // Size the FABs
246  //==========================================================
247  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
248  int ivar_idx = var_map[ivar];
249  Box domain = geom.Domain();
250  auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
251  domain.convert(ixtype);
252 
253  // NOTE: Ghost cells needed for idx type mismatch between mask and data (do_upwind)
254  IntVect ng_vect(0);
255  //IntVect ng_vect(1,1,0);
256  Box gdom(domain); gdom.grow(ng_vect);
257  Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
258  realbdy_interior_bxs_xy(gdom, domain, width,
259  bx_xlo, bx_xhi,
260  bx_ylo, bx_yhi,
261  ng_vect, true);
262 
263  // Size the FABs
264  if (ivar == ivarU) {
265  U_xlo.resize(bx_xlo,1,The_Async_Arena()); U_xhi.resize(bx_xhi,1,The_Async_Arena());
266  U_ylo.resize(bx_ylo,1,The_Async_Arena()); U_yhi.resize(bx_yhi,1,The_Async_Arena());
267  } else if (ivar == ivarV) {
268  V_xlo.resize(bx_xlo,1,The_Async_Arena()); V_xhi.resize(bx_xhi,1,The_Async_Arena());
269  V_ylo.resize(bx_ylo,1,The_Async_Arena()); V_yhi.resize(bx_yhi,1,The_Async_Arena());
270  } else if (ivar == ivarT){
271  T_xlo.resize(bx_xlo,1,The_Async_Arena()); T_xhi.resize(bx_xhi,1,The_Async_Arena());
272  T_ylo.resize(bx_ylo,1,The_Async_Arena()); T_yhi.resize(bx_yhi,1,The_Async_Arena());
273  } else {
274  continue;
275  }
276  } // ivar
277 
278 
279  // NOTE: These operations use the BDY FABS and RHO. The
280  // use of RHO to go from PRIM -> CONS requires that
281  // these operations be LOCAL. So we have allocated
282  // enough space to do global operations (1 rank) but
283  // will fill a subset of that data that the rank owns.
284 
285  // Populate FABs from bdy interpolation (primitive vars)
286  //==========================================================
287  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
288  int ivar_idx = var_map[ivar];
289  Box domain = geom.Domain();
290  auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
291  domain.convert(ixtype);
292  const auto& dom_lo = lbound(domain);
293  const auto& dom_hi = ubound(domain);
294 
295  // BndryReg idx and limiting
296  int bdy_comp = bnd_map[ivar];
297  const auto& dom_cc_lo = lbound(geom.Domain());
298  const auto& dom_cc_hi = ubound(geom.Domain());
299 
300 #ifdef _OPENMP
301 #pragma omp parallel if (Gpu::notInLaunchRegion())
302 #endif
303  for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
304  // NOTE: Ghost cells needed for idx type mismatch between mask and data (do_upwind)
305  IntVect ng_vect(0);
306  //IntVect ng_vect(1,1,0);
307  Box gtbx = grow(mfi.tilebox(ixtype.toIntVect()),ng_vect);
308  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
309  realbdy_interior_bxs_xy(gtbx, domain, width,
310  tbx_xlo, tbx_xhi,
311  tbx_ylo, tbx_yhi,
312  ng_vect, true);
313 
314  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
315  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
316  if (ivar == ivarU) {
317  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
318  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
319  } else if (ivar == ivarV) {
320  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
321  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
322  } else if (ivar == ivarT){
323  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
324  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
325  } else {
326  continue;
327  }
328 
329  // Boundary data at fixed time intervals
330  const auto& bdatxlo_n = bdy_data_xlo[n_time ][ivar].const_array();
331  const auto& bdatxlo_np1 = bdy_data_xlo[n_time_p1][ivar].const_array();
332  const auto& bdatxhi_n = bdy_data_xhi[n_time ][ivar].const_array();
333  const auto& bdatxhi_np1 = bdy_data_xhi[n_time_p1][ivar].const_array();
334  const auto& bdatylo_n = bdy_data_ylo[n_time ][ivar].const_array();
335  const auto& bdatylo_np1 = bdy_data_ylo[n_time_p1][ivar].const_array();
336  const auto& bdatyhi_n = bdy_data_yhi[n_time ][ivar].const_array();
337  const auto& bdatyhi_np1 = bdy_data_yhi[n_time_p1][ivar].const_array();
338 
339  // Current density to convert to conserved vars
340  Array4<Real> r_arr = S_cur_data[IntVars::cons].array(mfi);
341 
342  // Limiting offset
343  int offset = width - 1;
344 
345  // Populate with interpolation (protect from ghost cells)
346  ParallelFor(tbx_xlo, tbx_xhi,
347  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
348  {
349  int ii = std::max(i , dom_lo.x); ii = std::min(ii, dom_lo.x+offset);
350  int jj = std::max(j , dom_lo.y); jj = std::min(jj, dom_hi.y);
351 
352  Real rho_interp;
353  if (ivar==ivarU) {
354  rho_interp = myhalf * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
355  } else if (ivar==ivarV) {
356  rho_interp = myhalf * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
357  } else {
358  rho_interp = r_arr(i,j,k);
359  }
360 
361  if (bdatxlo) {
362  int ii2 = std::min(std::max(i , dom_cc_lo.x), dom_cc_hi.x);
363  int jj2 = std::min(std::max(j , dom_cc_lo.y), dom_cc_hi.y);
364  arr_xlo(i,j,k) = rho_interp * bdatxlo(ii2,jj2,k,bdy_comp);
365  } else {
366  arr_xlo(i,j,k) = rho_interp * ( oma * bdatxlo_n (ii,jj,k,0)
367  + alpha * bdatxlo_np1(ii,jj,k,0) );
368  }
369  },
370  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
371  {
372  int ii = std::max(i , dom_hi.x-offset); ii = std::min(ii, dom_hi.x);
373  int jj = std::max(j , dom_lo.y); jj = std::min(jj, dom_hi.y);
374 
375  Real rho_interp;
376  if (ivar==ivarU) {
377  rho_interp = myhalf * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
378  } else if (ivar==ivarV) {
379  rho_interp = myhalf * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
380  } else {
381  rho_interp = r_arr(i,j,k);
382  }
383 
384  if (bdatxhi) {
385  int ii2 = std::min(std::max(i , dom_cc_lo.x), dom_cc_hi.x);
386  int jj2 = std::min(std::max(j , dom_cc_lo.y), dom_cc_hi.y);
387  arr_xhi(i,j,k) = rho_interp * bdatxhi(ii2,jj2,k,bdy_comp);
388  } else {
389  arr_xhi(i,j,k) = rho_interp * ( oma * bdatxhi_n (ii,jj,k,0)
390  + alpha * bdatxhi_np1(ii,jj,k,0) );
391  }
392  });
393 
394  ParallelFor(tbx_ylo, tbx_yhi,
395  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
396  {
397  int ii = std::max(i , dom_lo.x); ii = std::min(ii, dom_hi.x);
398  int jj = std::max(j , dom_lo.y); jj = std::min(jj, dom_lo.y+offset);
399 
400  Real rho_interp;
401  if (ivar==ivarU) {
402  rho_interp = myhalf * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
403  } else if (ivar==ivarV) {
404  rho_interp = myhalf * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
405  } else {
406  rho_interp = r_arr(i,j,k);
407  }
408 
409  if (bdatylo) {
410  int ii2 = std::min(std::max(i , dom_cc_lo.x), dom_cc_hi.x);
411  int jj2 = std::min(std::max(j , dom_cc_lo.y), dom_cc_hi.y);
412  arr_ylo(i,j,k) = rho_interp * bdatylo(ii2,jj2,k,bdy_comp);
413  } else {
414  arr_ylo(i,j,k) = rho_interp * ( oma * bdatylo_n (ii,jj,k,0)
415  + alpha * bdatylo_np1(ii,jj,k,0) );
416  }
417  },
418  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
419  {
420  int ii = std::max(i , dom_lo.x); ii = std::min(ii, dom_hi.x);
421  int jj = std::max(j , dom_hi.y-offset); jj = std::min(jj, dom_hi.y);
422 
423  Real rho_interp;
424  if (ivar==ivarU) {
425  rho_interp = myhalf * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
426  } else if (ivar==ivarV) {
427  rho_interp = myhalf * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
428  } else {
429  rho_interp = r_arr(i,j,k);
430  }
431 
432  if (bdatyhi) {
433  int ii2 = std::min(std::max(i , dom_cc_lo.x), dom_cc_hi.x);
434  int jj2 = std::min(std::max(j , dom_cc_lo.y), dom_cc_hi.y);
435  arr_yhi(i,j,k) = rho_interp * bdatyhi(ii2,jj2,k,bdy_comp);
436  } else {
437  arr_yhi(i,j,k) = rho_interp * ( oma * bdatyhi_n (ii,jj,k,0)
438  + alpha * bdatyhi_np1(ii,jj,k,0) );
439  }
440  });
441  } // mfi
442  } // ivar
443 
444 
445  // Compute RHS in relaxation region
446  //==========================================================
447  auto dx = geom.CellSizeArray();
448  auto ProbLo = geom.ProbLoArray();
449  auto ProbHi = geom.ProbHiArray();
450 
451  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
452  int ivar_idx = ivar_map[ivar];
453  int icomp = comp_map[ivar];
454 
455  Box domain = geom.Domain();
456  domain.convert(S_cur_data[ivar_idx].boxArray().ixType());
457  IntVect ng_vect(0);
458 
459 #ifdef _OPENMP
460 #pragma omp parallel if (Gpu::notInLaunchRegion())
461 #endif
462  for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi)
463  {
464  Box tbx = mfi.tilebox();
465  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
466  realbdy_interior_bxs_xy(tbx, domain, width,
467  tbx_xlo, tbx_xhi,
468  tbx_ylo, tbx_yhi,
469  ng_vect);
470 
471  Array4<Real> rhs_arr; Array4<Real> data_arr;
472  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
473  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
474  if (ivar == ivarU) {
475  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
476  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
477  rhs_arr = S_rhs[IntVars::xmom].array(mfi);
478  data_arr = S_cur_data[IntVars::xmom].array(mfi);
479  } else if (ivar == ivarV) {
480  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
481  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
482  rhs_arr = S_rhs[IntVars::ymom].array(mfi);
483  data_arr = S_cur_data[IntVars::ymom].array(mfi);
484  } else if (ivar == ivarT){
485  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
486  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
487  rhs_arr = S_rhs[IntVars::cons].array(mfi);
488  data_arr = S_cur_data[IntVars::cons].array(mfi);
489  } else {
490  continue;
491  }
492 
494  width, dx, ProbLo, ProbHi, F1,
495  tbx_xlo , tbx_xhi , tbx_ylo , tbx_yhi ,
496  arr_xlo , arr_xhi , arr_ylo , arr_yhi ,
497  data_arr, rhs_arr);
498  } // mfi
499  } // ivar
500 
501  // Set normal velocity RHS at the boundary
502  //==========================================================
503  Box domain = geom.Domain();
504  Box domainx = convert(domain, IntVect(1,0,0));
505  Box domainy = convert(domain, IntVect(0,1,0));
506 
507  int ilo = domainx.smallEnd(0);
508  int ihi = domainx.bigEnd(0);
509  int jlo = domainy.smallEnd(1);
510  int jhi = domainy.bigEnd(1);
511 
512 #ifdef _OPENMP
513 #pragma omp parallel if (Gpu::notInLaunchRegion())
514 #endif
515  for (MFIter mfi(S_cur_data[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
516  {
517  Box tbx = mfi.nodaltilebox(0);
518  Box tbx_lo, tbx_hi;
519  if (tbx.smallEnd(0) == ilo) {
520  tbx_lo = makeSlab(tbx,0,ilo);
521  }
522  if (tbx.bigEnd(0) == ihi) {
523  tbx_hi = makeSlab(tbx,0,ihi);
524  }
525 
526  Box tby = mfi.nodaltilebox(1);
527  Box tby_lo, tby_hi;
528  if (tby.smallEnd(1) == jlo) {
529  tby_lo = makeSlab(tby,1,jlo);
530  }
531  if (tby.bigEnd(1) == jhi) {
532  tby_hi = makeSlab(tby,1,jhi);
533  }
534 
535  Array4<Real> rhs_xmom = S_rhs[IntVars::xmom].array(mfi);
536  Array4<Real> rhs_ymom = S_rhs[IntVars::ymom].array(mfi);
537 
538  Array4<const Real> rhs_cons = S_rhs[IntVars::cons].const_array(mfi);
539  Array4<const Real> cons_arr = S_cur_data[IntVars::cons].const_array(mfi);
540 
541  const auto& bdatxlo_n = bdy_data_xlo[n_time ][ivarU].const_array();
542  const auto& bdatxlo_np1 = bdy_data_xlo[n_time_p1][ivarU].const_array();
543  const auto& bdatxhi_n = bdy_data_xhi[n_time ][ivarU].const_array();
544  const auto& bdatxhi_np1 = bdy_data_xhi[n_time_p1][ivarU].const_array();
545 
546  const auto& bdatylo_n = bdy_data_ylo[n_time ][ivarV].const_array();
547  const auto& bdatylo_np1 = bdy_data_ylo[n_time_p1][ivarV].const_array();
548  const auto& bdatyhi_n = bdy_data_yhi[n_time ][ivarV].const_array();
549  const auto& bdatyhi_np1 = bdy_data_yhi[n_time_p1][ivarV].const_array();
550 
551  ParallelFor(tbx_lo, tbx_hi,
552  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
553  {
554  Real rho_tend = rhs_cons(i,j,k);
555  Real rho_val = Real(0.5) * (cons_arr(i,j,k) + cons_arr(i-1,j,k));
556  Real u_tend, u_val;
557  if (btenxlo) {
558  u_tend = btenxlo(i,j,k,BCVars::xvel_bc);
559  u_val = bdatxlo(i,j,k,BCVars::xvel_bc);
560  } else {
561  u_tend = (bdatxlo_np1(i,j,k) - bdatxlo_n(i,j,k)) / bdy_time_interval;
562  u_val = oma * bdatxlo_n(i,j,k) + alpha * bdatxlo_np1(i,j,k);
563  }
564  rhs_xmom(i,j,k) = rho_val * u_tend + u_val * rho_tend;
565  },
566  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
567  {
568  Real rho_tend = rhs_cons(i,j,k);
569  Real rho_val = Real(0.5) * (cons_arr(i,j,k) + cons_arr(i-1,j,k));
570  Real u_tend, u_val;
571  if (btenxhi) {
572  u_tend = btenxhi(i,j,k,BCVars::xvel_bc);
573  u_val = bdatxhi(i,j,k,BCVars::xvel_bc);
574  } else {
575  u_tend = (bdatxhi_np1(i,j,k) - bdatxhi_n(i,j,k)) / bdy_time_interval;
576  u_val = oma * bdatxhi_n(i,j,k) + alpha * bdatxhi_np1(i,j,k);
577  }
578  rhs_xmom(i,j,k) = rho_val * u_tend + u_val * rho_tend;
579  });
580 
581  ParallelFor(tby_lo, tby_hi,
582  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
583  {
584  Real rho_tend = rhs_cons(i,j,k);
585  Real rho_val = Real(0.5) * (cons_arr(i,j,k) + cons_arr(i,j-1,k));
586  Real v_tend, v_val;
587  if (btenylo) {
588  v_tend = btenylo(i,j,k,BCVars::yvel_bc);
589  v_val = bdatylo(i,j,k,BCVars::yvel_bc);
590  } else {
591  v_tend = (bdatylo_np1(i,j,k) - bdatylo_n(i,j,k)) / bdy_time_interval;
592  v_val = oma * bdatylo_n(i,j,k) + alpha * bdatylo_np1(i,j,k);
593  }
594  rhs_ymom(i,j,k) = rho_val * v_tend + v_val * rho_tend;
595  },
596  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
597  {
598  Real rho_tend = rhs_cons(i,j,k);
599  Real rho_val = Real(0.5) * (cons_arr(i,j,k) + cons_arr(i,j-1,k));
600  Real v_tend, v_val;
601  if (btenyhi) {
602  v_tend = btenyhi(i,j,k,BCVars::yvel_bc);
603  v_val = bdatyhi(i,j,k,BCVars::yvel_bc);
604  } else {
605  v_tend = (bdatyhi_np1(i,j,k) - bdatyhi_n(i,j,k)) / bdy_time_interval;
606  v_val = oma * bdatyhi_n(i,j,k) + alpha * bdatyhi_np1(i,j,k);
607  }
608  rhs_ymom(i,j,k) = rho_val * v_tend + v_val * rho_tend;
609  });
610  } // mfi
611 }
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
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
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 &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:225
@ RhoTheta_bc_comp
Definition: ERF_IndexDefines.H:88
@ yvel_bc
Definition: ERF_IndexDefines.H:103
@ xvel_bc
Definition: ERF_IndexDefines.H:102
@ U
Definition: ERF_IndexDefines.H:123
@ NumTypes
Definition: ERF_IndexDefines.H:127
@ T
Definition: ERF_IndexDefines.H:125
@ V
Definition: ERF_IndexDefines.H:124
@ xvel
Definition: ERF_IndexDefines.H:175
@ cons
Definition: ERF_IndexDefines.H:174
@ yvel
Definition: ERF_IndexDefines.H:176
real(kind=kind_phys), parameter, public alpha
Definition: ERF_module_mp_wsm6.F90:44
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 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
32 {
33  AMREX_ALWAYS_ASSERT(bx.ixType() == domain.ixType());
34 
35  //==================================================================
36  // NOTE: X-face boxes take ownership of the overlapping region.
37  // With exterior ghost cells (ng_vect != 0), the x-face
38  // boxes will have exterior ghost cells in both x & y.
39  //==================================================================
40 
41  // Domain bounds without ghost cells
42  const auto& dom_lo = lbound(domain);
43  const auto& dom_hi = ubound(domain);
44 
45  // Four boxes matching the domain
46  Box gdom_xlo(domain); Box gdom_xhi(domain);
47  Box gdom_ylo(domain); Box gdom_yhi(domain);
48 
49  // Trim the boxes to only include internal ghost cells
50  gdom_xlo.setBig(0,dom_lo.x+width-1); gdom_xhi.setSmall(0,dom_hi.x-width+1);
51  gdom_ylo.setBig(1,dom_lo.y+width-1); gdom_yhi.setSmall(1,dom_hi.y-width+1);
52 
53  // Remove overlapping corners from y-face boxes
54  gdom_ylo.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_ylo.setBig(0,gdom_xhi.smallEnd(0)-1);
55  gdom_yhi.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_yhi.setBig(0,gdom_xhi.smallEnd(0)-1);
56 
57  // Grow boxes to get external ghost cells only
58  gdom_xlo.growLo(0,ng_vect[0]); gdom_xhi.growHi(0,ng_vect[0]);
59  gdom_xlo.grow (1,ng_vect[1]); gdom_xhi.grow (1,ng_vect[1]);
60  gdom_ylo.growLo(1,ng_vect[1]); gdom_yhi.growHi(1,ng_vect[1]);
61 
62  // Grow boxes to get internal ghost cells
63  if (get_int_ng) {
64  gdom_xlo.growHi(0,ng_vect[0]); gdom_xhi.growLo(0,ng_vect[0]);
65  gdom_ylo.grow (0,ng_vect[0]); gdom_yhi.grow (0,ng_vect[0]);
66  gdom_ylo.growHi(1,ng_vect[1]); gdom_yhi.growLo(1,ng_vect[1]);
67  }
68 
69  // Populate everything
70  bx_xlo = (bx & gdom_xlo);
71  bx_xhi = (bx & gdom_xhi);
72  bx_ylo = (bx & gdom_ylo);
73  bx_yhi = (bx & gdom_yhi);
74 }

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