Skip to content

Commit

Permalink
New feature in EB: Add regular coarse levels
Browse files Browse the repository at this point in the history
This feature could be useful when the EB is so thin that it is only resolved
at the finest level and the usual way of generate coarse levels fails with
multi-cut or multi-value cells. When the EB is fully refined to the finest
level and there is no subcycling, no EB information is really needed at
coarse levels. Note that tagging cut cells for refinement no longer works if
that level has all regular EB.
  • Loading branch information
WeiqunZhang committed Jul 9, 2023
1 parent b41f85a commit 2f0b22e
Show file tree
Hide file tree
Showing 8 changed files with 77 additions and 0 deletions.
4 changes: 4 additions & 0 deletions Src/EB/AMReX_EB2.H
Original file line number Diff line number Diff line change
Expand Up @@ -56,6 +56,7 @@ public:
[[nodiscard]] virtual const Geometry& getGeometry (const Box& domain) const = 0;
[[nodiscard]] virtual const Box& coarsestDomain () const = 0;
virtual void addFineLevels (int num_new_fine_levels) = 0;
virtual void addRegularCoarseLevels (int num_new_coarse_levels) = 0;

protected:
static AMREX_EXPORT Vector<std::unique_ptr<IndexSpace> > m_instance;
Expand Down Expand Up @@ -88,6 +89,7 @@ public:
return m_geom.back().Domain();
}
void addFineLevels (int num_new_fine_levels) final;
void addRegularCoarseLevels (int num_new_coarse_levels) final;

using F = typename G::FunctionType;

Expand Down Expand Up @@ -147,6 +149,8 @@ int maxCoarseningLevel (IndexSpace const* ebis, const Geometry& geom);

void addFineLevels (int num_new_fine_levels);

void addRegularCoarseLevels (int num_new_coarse_levels);

}

#endif
8 changes: 8 additions & 0 deletions Src/EB/AMReX_EB2.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -240,6 +240,14 @@ void addFineLevels (int num_new_fine_levels)
}
}

void addRegularCoarseLevels (int num_new_coarse_levels)
{
auto *p = const_cast<IndexSpace*>(TopIndexSpace());
if (p) {
p->addRegularCoarseLevels(num_new_coarse_levels);
}
}

void
BuildFromChkptFile (std::string const& fname,
const Geometry& geom, int required_coarsening_level,
Expand Down
42 changes: 42 additions & 0 deletions Src/EB/AMReX_EB2_IndexSpaceI.H
Original file line number Diff line number Diff line change
Expand Up @@ -106,3 +106,45 @@ IndexSpaceImp<G>::addFineLevels (int num_new_fine_levels)
m_domain.insert(m_domain.begin(), fine_isp.m_domain.begin(), fine_isp.m_domain.end());
m_ngrow.insert(m_ngrow.begin(), fine_isp.m_ngrow.begin(), fine_isp.m_ngrow.end());
}

template <typename G>
void
IndexSpaceImp<G>::addRegularCoarseLevels (int num_new_coarse_levels)
{
if (num_new_coarse_levels <= 0) { return; }

auto nlevs_old = int(m_gslevel.size());
int nlevs_new = nlevs_old + num_new_coarse_levels;

Vector<GShopLevel<G>> new_gslevel;
new_gslevel.reserve(nlevs_new);

Vector<Geometry> new_geom;
new_geom.reserve(nlevs_new);

Vector<Box> new_domain;
new_domain.reserve(nlevs_new);

Vector<int> new_ngrow;
new_ngrow.reserve(nlevs_new);

for (int ilev = 0; ilev < num_new_coarse_levels; ++ilev) {
int rr = 1 << (num_new_coarse_levels-ilev); // 2^(num_new_coarse_levels-ilev)
new_geom.push_back(amrex::coarsen(m_geom[0], rr));
new_domain.push_back(new_geom.back().Domain());
new_ngrow.push_back(m_ngrow[0]);
new_gslevel.push_back(GShopLevel<G>::makeAllRegular(this, new_geom.back()));
}

for (int ilev = 0; ilev < nlevs_old; ++ilev) {
new_gslevel.emplace_back(std::move(m_gslevel[ilev]));
new_geom.push_back (m_geom [ilev]);
new_domain.push_back(m_domain[ilev]);
new_ngrow.push_back (m_ngrow [ilev]);
}

std::swap(new_gslevel, m_gslevel);
std::swap(new_geom , m_geom);
std::swap(new_domain , m_domain);
std::swap(new_ngrow , m_ngrow);
}
1 change: 1 addition & 0 deletions Src/EB/AMReX_EB2_IndexSpace_STL.H
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,7 @@ public:
return m_geom.back().Domain();
}
void addFineLevels (int num_new_fine_levels) final;
void addRegularCoarseLevels (int num_new_coarse_levels) final;

private:

Expand Down
6 changes: 6 additions & 0 deletions Src/EB/AMReX_EB2_IndexSpace_STL.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -88,4 +88,10 @@ IndexSpaceSTL::addFineLevels (int /*num_new_fine_levels*/)
amrex::Abort("IndexSpaceSTL::addFineLevels: todo");
}

void
IndexSpaceSTL::addRegularCoarseLevels (int /*num_new_coarse_levels*/)
{
amrex::Abort("IndexSpaceSTL::addRegularCoarseLevels: todo");
}

}
1 change: 1 addition & 0 deletions Src/EB/AMReX_EB2_IndexSpace_chkpt_file.H
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,7 @@ public:
return m_geom.back().Domain();
}
void addFineLevels (int num_new_fine_levels) final;
void addRegularCoarseLevels (int num_new_coarse_levels) final;

private:

Expand Down
6 changes: 6 additions & 0 deletions Src/EB/AMReX_EB2_IndexSpace_chkpt_file.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -83,4 +83,10 @@ IndexSpaceChkptFile::addFineLevels (int /*num_new_fine_levels*/)
amrex::Abort("IndexSpaceChkptFile::addFineLevels: not supported");
}

void
IndexSpaceChkptFile::addRegularCoarseLevels (int /*num_new_coarse_levels*/)
{
amrex::Abort("IndexSpaceChkptFile::addRegularCoarseLevels: not supported");
}

}
9 changes: 9 additions & 0 deletions Src/EB/AMReX_EB2_Level.H
Original file line number Diff line number Diff line change
Expand Up @@ -101,6 +101,15 @@ public:
GShopLevel (IndexSpace const* is, const Geometry& geom);
void define_fine (G const& gshop, const Geometry& geom,
int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt);

static GShopLevel<G>
makeAllRegular(IndexSpace const* is, const Geometry& geom)
{
GShopLevel<G> r(is, geom);
r.m_allregular = true;
r.m_ok = true;
return r;
}
};

template <typename G>
Expand Down

0 comments on commit 2f0b22e

Please sign in to comment.