Skip to content

Commit

Permalink
comments
Browse files Browse the repository at this point in the history
  • Loading branch information
RevathiJambunathan committed Feb 26, 2024
1 parent d72fd6e commit be9fe48
Showing 1 changed file with 13 additions and 5 deletions.
18 changes: 13 additions & 5 deletions Source/BoundaryConditions/PML.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -560,25 +560,33 @@ PML::PML (const int lev, const BoxArray& grid_ba, const DistributionMapping& gri
m_geom(geom),
m_cgeom(cgeom)
{
// When `do_pml_in_domain` is true, the PML overlap with the last `ncell` of the physical domain
// (instead of extending `ncell` outside of the physical domain)
// In order to implement this, a reduced domain is created here (decreased by ncells in all direction)
// and passed to `MakeBoxArray`, which surrounds it by PML boxes
// (thus creating the PML boxes at the right position, where they overlap with the original domain)
// When `do_pml_in_domain` is true, the PML overlap with the last `ncell` of the physical domain or fine patch(es)
// (instead of extending `ncell` outside of the physical domain or fine patch(es))
// In order to implement this, we define a new reduced Box Array ensuring that it does not
// include ncells from the edges of the physical domain or fine patch.
// (thus creating the PML boxes at the right position, where they overlap with the original domain or fine patch(es))

BoxArray grid_ba_reduced = grid_ba;
if (do_pml_in_domain) {
BoxList bl = grid_ba.boxList();
// Here we loop over all the boxes in the original grid_ba BoxArray
// For each box, we find if its in the edge (or boundary), and the size of those boxes are decreased by ncell
for (auto& b : bl) {
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
if (do_pml_Lo[idim]) {
// Get neighboring box on lower side in direction idim and check if it intersects with any of the boxes
// in grid_ba. If no intersection, then the box, b, in the boxlist, is in the edge and we decrase
// the size by ncells using growLo(idim,-ncell)
Box const& bb = amrex::adjCellLo(b, idim);
if ( ! grid_ba.intersects(bb) ) {
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(b.length(idim) > ncell, " box length must be greater that pml size");
b.growLo(idim, -ncell);
}
}
if (do_pml_Hi[idim]) {
// Get neighboring box on higher side in direction idim and check if it intersects with any of the boxes
// in grid_ba. If no intersection, then the box, b, in the boxlist, is in the edge and we decrase
// the size by ncells using growHi(idim,-ncell)
Box const& bb = amrex::adjCellHi(b, idim);
if ( ! grid_ba.intersects(bb) ) {
WARPX_ALWAYS_ASSERT_WITH_MESSAGE(b.length(idim) > ncell, " box length must be greater that pml size");
Expand Down

0 comments on commit be9fe48

Please sign in to comment.