Skip to content
Snippets Groups Projects
Commit c34f1855 authored by Praetorius, Simon's avatar Praetorius, Simon
Browse files

Merge branch 'issue/remove_marker' into 'master'

add functionality to remove an added marker

See merge request !46
parents 3fec45ce 2a0d6e7f
No related branches found
No related tags found
1 merge request!46add functionality to remove an added marker
......@@ -253,10 +253,11 @@ namespace AMDiS
public:
/// Implementation of \ref ProblemStatBase::solve
void solve(AdaptInfo& adaptInfo,
bool createMatrixData = true,
bool storeMatrixData = false) override;
/// Implementation of \ref StandardProblemIteration::oneIteration.
Flag oneIteration(AdaptInfo& adaptInfo, Flag toDo = FULL_ITERATION) override
{
return StandardProblemIteration::oneIteration(adaptInfo, toDo);
}
/// Implementation of \ref ProblemStatBase::buildAfterCoarse
void buildAfterAdapt(AdaptInfo& adaptInfo,
......@@ -271,6 +272,26 @@ namespace AMDiS
buildAfterAdapt(adaptInfo, Flag{0}, true, true);
}
/// Implementation of \ref ProblemStatBase::solve
void solve(AdaptInfo& adaptInfo,
bool createMatrixData = true,
bool storeMatrixData = false) override;
/// Implementation of \ref ProblemStatBase::estimate.
void estimate(AdaptInfo& adaptInfo) override { /* do nothing. */ }
/// Implementation of \ref ProblemStatBase::refineMesh.
Flag adaptGrid(AdaptInfo& adaptInfo) override;
/// Implementation of \ref ProblemStatBase::markElements.
Flag markElements(AdaptInfo& adaptInfo) override;
/// Uniform global grid coarsening by up to n level
Flag globalCoarsen(int n) override;
/// Uniform global refinement by n level
Flag globalRefine(int n) override;
/// Writes output files.
void writeFiles(AdaptInfo& adaptInfo, bool force = false);
......@@ -335,14 +356,18 @@ namespace AMDiS
public: // set-methods
/// Set a new linear solver for the problem
void setSolver(std::shared_ptr<LinearSolverType> const& solver)
void setSolver(std::shared_ptr<LinearSolverType> solver)
{
linearSolver_ = solver;
}
void setSolver(LinearSolverType& solver)
/// Wrap the solver reference into a non-destroying shared_ptr, or move it
/// into a new shared_ptr if it is a temporary.
template <class S,
REQUIRES(std::is_base_of<LinearSolverType, remove_cvref_t<S>>::value)>
void setSolver(S&& solver)
{
setSolver(Dune::stackobject_to_shared_ptr(solver));
setSolver(Dune::wrap_or_move(FWD(solver)));
}
......@@ -357,21 +382,44 @@ namespace AMDiS
createFileWriter();
}
/// Wrap the grid into a non-destroying shared_ptr
void setGrid(Grid& grid)
{
setGrid(Dune::stackobject_to_shared_ptr(grid));
}
void addMarker(std::shared_ptr<Marker<Grid>> const& marker)
/// Store the shared_ptr and the name of the marker in the problem
/**
* Note: multiple markers can be added but must have different names
**/
void addMarker(std::shared_ptr<Marker<Grid>> marker)
{
marker_.push_back(marker);
auto it = marker_.emplace(marker->name(), std::move(marker));
if (marker_.size() > 1)
marker_.back()->setMaximumMarking(true);
it.first->second->setMaximumMarking(true);
}
void addMarker(Marker<Grid>& marker)
/// Wrap the reference into a non-destroying shared_ptr or move it into a
/// new shared_ptr if it is a temporary.
template <class M,
REQUIRES(std::is_base_of<Marker<Grid>, remove_cvref_t<M>>::value)>
void addMarker(M&& marker)
{
addMarker(Dune::stackobject_to_shared_ptr(marker));
addMarker(Dune::wrap_or_move(FWD(marker)));
}
/// Remove a marker with the given name from the problem
void removeMarker(std::string name)
{
std::size_t num = marker_.erase(name);
test_warning(num == 1, "A marker with the given name '{}' does not exist.", name);
}
/// Remove a marker from the problem
void removeMarker(Marker<Grid> const& marker)
{
removeMarker(marker.name());
}
......@@ -390,11 +438,6 @@ namespace AMDiS
initGlobalBasis();
}
void adoptGrid(std::shared_ptr<Grid> const& grid)
{
adoptGrid(grid, std::make_shared<BoundaryManager<Grid>>(grid));
}
void adoptGrid(std::shared_ptr<Grid> const& grid,
std::shared_ptr<BoundaryManager<Grid>> const& boundaryManager)
{
......@@ -403,6 +446,11 @@ namespace AMDiS
Parameters::get(name_ + "->mesh", gridName_);
}
void adoptGrid(std::shared_ptr<Grid> const& grid)
{
adoptGrid(grid, std::make_shared<BoundaryManager<Grid>>(grid));
}
private:
void createGlobalBasisImpl(std::true_type);
......@@ -410,30 +458,6 @@ namespace AMDiS
void initGlobalBasis();
public: // implementation of iteration interface methods
/// Implementation of \ref StandardProblemIteration::oneIteration.
Flag oneIteration(AdaptInfo& adaptInfo, Flag toDo = FULL_ITERATION) override
{
return StandardProblemIteration::oneIteration(adaptInfo, toDo);
}
/// Implementation of \ref ProblemStatBase::estimate.
void estimate(AdaptInfo& adaptInfo) override { /* do nothing. */ }
/// Implementation of \ref ProblemStatBase::refineMesh.
Flag adaptGrid(AdaptInfo& adaptInfo) override;
/// Implementation of \ref ProblemStatBase::markElements.
Flag markElements(AdaptInfo& adaptInfo) override;
/// Uniform global grid coarsening by up to n level
Flag globalCoarsen(int n) override;
/// Uniform global refinement by n level
Flag globalRefine(int n) override;
private:
/// Name of this problem.
std::string name_;
......@@ -454,7 +478,7 @@ namespace AMDiS
std::list<std::shared_ptr<FileWriterInterface>> filewriter_;
/// Pointer to the adaptation markers
std::list<std::shared_ptr<Marker<Grid>>> marker_;
std::map<std::string, std::shared_ptr<Marker<Grid>>> marker_;
/// Pointer to the estimators for this problem
// std::vector<Estimator*> estimator;
......
......@@ -217,14 +217,8 @@ void ProblemStat<Traits>::createMarker()
std::string tp = to_string(treePath);
auto newMarker
= EstimatorMarker<Grid>::createMarker(componentName, tp, estimates_[tp], grid_);
if (newMarker)
marker_.push_back(std::move(newMarker));
// If there is more than one marker, and all components are defined
// on the same grid, the maximum marking has to be enabled.
// TODO: needs to be checked!
if (marker_.size() > 1)
marker_.back()->setMaximumMarking(true);
assert(bool(newMarker));
this->addMarker(std::move(newMarker));
});
}
......@@ -348,7 +342,7 @@ markElements(AdaptInfo& adaptInfo)
Flag markFlag = 0;
for (auto& currentMarker : marker_)
markFlag |= currentMarker->markGrid(adaptInfo);
markFlag |= currentMarker.second->markGrid(adaptInfo);
msg("markElements needed {} seconds", t.elapsed());
......
......@@ -53,5 +53,10 @@ int main(int argc, char** argv)
prob.solution().interpolate(markerFunc);
prob.writeFiles(adaptInfo);
prob.removeMarker(marker.name());
prob.addMarker(marker);
prob.removeMarker(marker);
return report_errors();
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment