Commit 681e70ef authored by Praetorius, Simon's avatar Praetorius, Simon
Browse files

Constraints improvements

parent 34c64903
#pragma once
#include <list>
#include <memory>
#include <type_traits>
#include <amdis/typetree/TreeContainer.hpp>
#include <dune/functions/common/typeerasure.hh>
namespace AMDiS
#include <amdis/common/ConceptsBase.hpp>
#include <amdis/common/TypeTraits.hpp>
namespace AMDiS {
namespace Impl {
template <class Mat, class Sol, class Rhs>
struct BoundaryConditionDefinition
{
/// \brief Interface class for boundary conditions
/**
* Stores a boundary subset related to the boundary condition. See \ref BoundarySubset.
* \relates DirichletBoundaryCondition
* \relates PeriodicBoundaryCondition
* \tparam Mat Matrix
* \tparam Sol Vector of solution
* \tparam Rhs Vector of rhs
**/
template <class Mat, class Sol, class Rhs>
class BoundaryCondition
// Definition of the interface a BoundaryCondition must fulfill
struct Interface
{
public:
/// Default constructor.
BoundaryCondition() = default;
/// \brief Initialize the boundary condition.
/**
* This performs setup before starting the matrix and vector assembly independantly of the
* matrix, rhs or solution.
**/
virtual ~Interface() = default;
virtual void init() = 0;
virtual void apply(Mat&, Sol&, Rhs&) = 0;
};
/// \brief Apply the boundary condition to matrix and vector.
/**
* This is called after the matrix `A` and rhs `b` are assembled.
* Implementations may alter values or change the structure of `A`, `x` or `b`.
**/
virtual void apply(Mat& A, Sol& x, Rhs& b) = 0;
// Templatized implementation of the interface
template <class Impl>
struct Model : public Impl
{
using Impl::Impl;
void init() final { this->get().init(); }
void apply(Mat& A, Sol& x, Rhs& b) final { this->get().apply(A,x,b); }
};
template <class Mat, class Sol, class Rhs>
struct BCData
// The Concept definition of a BoundaryCondition
struct Concept
{
template <class RN, class CN>
using type = std::list<std::shared_ptr<BoundaryCondition<Mat, Sol, Rhs>>>;
template <class BC>
auto require(BC&& bc) -> decltype
(
bc.init(),
bc.apply(std::declval<Mat&>(), std::declval<Sol&>(), std::declval<Rhs&>())
);
};
/// Container for storing shared pointers to boundary conditions indexed by their position in
/// the basis tree, see \ref MatrixData.
template <class Mat, class Sol, class Rhs, class RB, class CB>
using BoundaryConditions
= TreeMatrix<BCData<Mat, Sol, Rhs>::template type,
typename RB::LocalView::Tree, typename CB::LocalView::Tree>;
using Base = Dune::Functions::TypeErasureBase<Interface, Model>;
};
} // end namespace Impl
/// \brief Interface class for boundary conditions
/**
* Stores a boundary subset related to the boundary condition. See \ref BoundarySubset.
* \relates DirichletBoundaryCondition
* \relates PeriodicBoundaryCondition
* \tparam Mat Matrix
* \tparam Sol Vector of solution
* \tparam Rhs Vector of rhs
**/
template <class Mat, class Sol, class Rhs>
class BoundaryCondition
: public Impl::BoundaryConditionDefinition<Mat,Sol,Rhs>::Base
{
using Definition = Impl::BoundaryConditionDefinition<Mat,Sol,Rhs>;
using Super = typename Definition::Base;
public:
/// \brief Constructor. Pass any type supporting the \ref BoundaryConditionInterface
template <class Impl, Dune::disableCopyMove<BoundaryCondition,Impl> = 0>
BoundaryCondition(Impl&& impl)
: Super{FWD(impl)}
{
static_assert(Concepts::models<typename Definition::Concept(Impl)>,
"Implementation does not model the BoundaryCondition concept.");
}
/// \brief Default Constructor.
BoundaryCondition() = default;
/// \brief Initialize the boundary condition.
/**
* This performs setup before starting the matrix and vector assembly independantly of the
* matrix, rhs or solution.
**/
void init()
{
this->asInterface().init();
}
/// \brief Apply the boundary condition to matrix and vector.
/**
* This is called after the matrix `A` and rhs `b` are assembled.
* Implementations may alter values or change the structure of `A`, `x` or `b`.
**/
void apply(Mat& A, Sol& x, Rhs& b)
{
this->asInterface().apply(A,x,b);
}
};
} // end namespace AMDiS
......@@ -32,22 +32,14 @@ namespace AMDiS
* application of boundary conditions can be symmetric if the matrix does
* support this symmetric modification.
*
* \tparam Mat Matrix
* \tparam Sol Vector of solution
* \tparam Rhs Vector of rhs
* \tparam RowSubBasis SubspaceBasis of the row FE space
* \tparam ColSubBasis SubspaceBasis of the column FE space
* \tparam SubBasis SubspaceBasis of the solution FE space
* \tparam ValueGridFct Type of the GridFunction representing the Dirichlet values.
**/
template <class Mat, class Sol, class Rhs, class RowSubBasis, class ColSubBasis,
class ValueGridFct>
template <class SubBasis, class ValueGridFct>
class DirichletBC
: public BoundaryCondition<Mat, Sol, Rhs>
{
using Super = BoundaryCondition<Mat, Sol, Rhs>;
// We assume CB's GridView to be the same as RB's
using GridView = typename RowSubBasis::GridView;
using GridView = typename SubBasis::GridView;
using Intersection = typename GridView::Intersection;
using Domain = typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
......@@ -57,78 +49,72 @@ namespace AMDiS
/// Constructor accepting subspacebases.
template <class Values,
REQUIRES(Concepts::Functor<Values, Range(Domain)>)>
DirichletBC(RowSubBasis rowBasis, ColSubBasis colBasis,
BoundarySubset<Intersection> boundarySubset, Values&& values)
: rowBasis_(std::move(rowBasis))
, colBasis_(std::move(colBasis))
DirichletBC(SubBasis basis, BoundarySubset<Intersection> boundarySubset, Values&& values)
: basis_(std::move(basis))
, boundarySubset_(std::move(boundarySubset))
, values_(FWD(values))
{}
/// Make a DirichletBC from a basis with treepath arguments
template <class Basis, class TP, class Values,
REQUIRES(Concepts::GlobalBasis<Basis>)>
DirichletBC(Basis const& basis, TP const& treePath,
BoundarySubset<Intersection> boundarySubset, Values&& values)
: DirichletBC(Dune::Functions::subspaceBasis(basis, treePath),
std::move(boundarySubset), FWD(values))
{}
/// Make a DirichletBC from a global basis
template <class Basis, class Values,
REQUIRES(Concepts::GlobalBasis<Basis>)>
DirichletBC(Basis const& basis, BoundarySubset<Intersection> boundarySubset, Values&& values)
: DirichletBC(Dune::Functions::subspaceBasis(basis, makeTreePath()),
std::move(boundarySubset), FWD(values))
{}
/// Fill \ref dirichletNodes_ with 1 or 0 whether DOF corresponds to the
/// boundary or not.
/**
* \see BoundaryCondition::init
**/
void init() override;
void init();
/// \brief Apply dirichlet BC to matrix and vector
/**
* Add a unit-row to the matrix and optionally delete the corresponding matrix-column.
* Uses a backend-specific implementation.
**/
void apply(Mat& matrix, Sol& solution, Rhs& rhs) override;
template <class Mat, class Sol, class Rhs>
void apply(Mat& matrix, Sol& solution, Rhs& rhs);
private:
RowSubBasis rowBasis_;
ColSubBasis colBasis_;
SubBasis basis_;
BoundarySubset<Intersection> boundarySubset_;
ValueGridFct values_;
std::vector<bool> dirichletNodes_;
};
/// Make a DirichletBC from subspacebases
template <class Mat, class Sol, class Rhs, class RB, class RTP, class CB, class CTP, class Values>
auto makeDirichletBC(Dune::Functions::SubspaceBasis<RB, RTP> rowBasis,
Dune::Functions::SubspaceBasis<CB, CTP> colBasis,
BoundarySubset<typename RB::GridView::Intersection> boundarySubset,
Values&& values)
{
using RowSubBasis = Dune::Functions::SubspaceBasis<RB, RTP>;
using ColSubBasis = Dune::Functions::SubspaceBasis<CB, CTP>;
using BcType = DirichletBC<Mat, Sol, Rhs, RowSubBasis, ColSubBasis, TYPEOF(values)>;
return BcType(std::move(rowBasis), std::move(colBasis), std::move(boundarySubset), FWD(values));
}
/// Make a DirichletBC from a row- and colbasis with treepath arguments
template <class Mat, class Sol, class Rhs, class RB, class RTP, class CB, class CTP, class Values,
REQUIRES(Concepts::GlobalBasis<RB>),
REQUIRES(Concepts::GlobalBasis<CB>)>
auto makeDirichletBC(RB const& rowBasis, RTP const& rowTreePath,
CB const& colBasis, CTP const& colTreePath,
BoundarySubset<typename RB::GridView::Intersection> boundarySubset,
Values&& values)
{
return makeDirichletBC<Mat,Sol,Rhs>(
Dune::Functions::subspaceBasis(rowBasis, rowTreePath),
Dune::Functions::subspaceBasis(colBasis, colTreePath),
std::move(boundarySubset), FWD(values));
}
/// Make a DirichletBC from a global row- and colbasis
template <class Mat, class Sol, class Rhs, class RB, class CB, class Values,
REQUIRES(Concepts::GlobalBasis<RB>),
REQUIRES(Concepts::GlobalBasis<CB>)>
auto makeDirichletBC(RB const& rowBasis, CB const& colBasis,
BoundarySubset<typename RB::GridView::Intersection> boundarySubset,
Values&& values)
{
return makeDirichletBC<Mat, Sol, Rhs>(
rowBasis, makeTreePath(), colBasis, makeTreePath(),
std::move(boundarySubset), FWD(values));
}
// deduction guides
// Make a DirichletBC from subspacebases
template <class B, class TP, class Values>
DirichletBC(Dune::Functions::SubspaceBasis<B, TP>,
BoundarySubset<typename B::GridView::Intersection>, Values const&)
-> DirichletBC<Dune::Functions::SubspaceBasis<B,TP>, Underlying_t<Values>>;
// Make a DirichletBC from a basis with treepath arguments
template <class B, class TP, class Values,
REQUIRES(Concepts::GlobalBasis<B>)>
DirichletBC(B const&, TP const&, BoundarySubset<typename B::GridView::Intersection>, Values const&)
-> DirichletBC<Dune::Functions::SubspaceBasis<B, TP>, Underlying_t<Values>>;
// Make a DirichletBC from a global basis
template <class B, class Values,
REQUIRES(Concepts::GlobalBasis<B>)>
DirichletBC(B const&, BoundarySubset<typename B::GridView::Intersection>, Values const&)
-> DirichletBC<Dune::Functions::SubspaceBasis<B,TYPEOF(makeTreePath())>, Underlying_t<Values>>;
} // end namespace AMDiS
......
......@@ -8,50 +8,36 @@
#include <amdis/LinearAlgebra.hpp>
namespace AMDiS {
namespace Impl {
/// Traverse all interior boundary DOFs and apply the functor f.
/// The Functor must have the signature
/// `void(int, typename Basis::LocalView, typename Basis::GridView::Intersection)`
template <class Basis, class F>
void forEachInteriorBoundaryDOF(Basis const& basis, F&& f)
template <class B, class V>
void DirichletBC<B, V>::
init()
{
auto localView = basis.localView();
auto seDOFs = Dune::Functions::subEntityDOFs(basis);
auto const& gridView = basis.gridView();
for (auto const& element : elements(gridView, typename BackendTraits<Basis>::PartitionSet{})) {
dirichletNodes_.assign(basis_.dimension(), false);
auto localView = basis_.localView();
auto seDOFs = Dune::Functions::subEntityDOFs(basis_);
auto const& gridView = basis_.gridView();
for (auto const& element : elements(gridView, typename BackendTraits<B>::PartitionSet{})) {
if (element.hasBoundaryIntersections()) {
localView.bind(element);
for(auto const& intersection: intersections(gridView, element))
if (intersection.boundary())
for(auto localIndex: seDOFs.bind(localView,intersection))
f(localIndex, localView, intersection);
dirichletNodes_[localView.index(localIndex)]
= dirichletNodes_[localView.index(localIndex)] || boundarySubset_(intersection);
}
}
}
} // end namespace Impl
template <class Mat, class Sol, class Rhs, class RB, class CB, class V>
void DirichletBC<Mat, Sol, Rhs, RB, CB, V>::
init()
{
using LV = typename CB::LocalView;
dirichletNodes_.assign(colBasis_.dimension(), false);
Impl::forEachInteriorBoundaryDOF(colBasis_,
[&](int localIndex, LV const& localView, Intersection const& intersection) {
dirichletNodes_[localView.index(localIndex)] = dirichletNodes_[localView.index(localIndex)]
|| boundarySubset_(intersection);
});
}
template <class Mat, class Sol, class Rhs, class RB, class CB, class V>
void DirichletBC<Mat, Sol, Rhs, RB, CB, V>::
template <class B, class V>
template <class Mat, class Sol, class Rhs>
void DirichletBC<B, V>::
apply(Mat& matrix, Sol& solution, Rhs& rhs)
{
AMDiS::interpolate(colBasis_, solution, values_, tag::defaulted{}, tag::defaulted{}, dirichletNodes_);
AMDiS::interpolate(basis_, solution, values_, tag::defaulted{}, tag::defaulted{}, dirichletNodes_);
dirichletBC(matrix, solution, rhs, dirichletNodes_);
}
......
......@@ -57,11 +57,9 @@ namespace AMDiS
* \tparam Basis RowBasis, must be the same as ColBasis
* \tparam TP Treepath to the node this constraint applies to
**/
template <class Mat, class Sol, class Rhs, class Basis, class TP>
template <class Basis, class TP>
class PeriodicBC
: public BoundaryCondition<Mat, Sol, Rhs>
{
using Super = BoundaryCondition<Mat, Sol, Rhs>;
using SubspaceBasis = Dune::Functions::SubspaceBasis<Basis, TP>;
public:
......@@ -85,7 +83,7 @@ namespace AMDiS
* true for all DOFs on the periodic boundary. Furthermode associations() is initialized and
* contains the index of the DOF associated to a DOF on the periodic boundary.
**/
void init() override;
void init();
/// Move matrix rows (and columns) of dependant DOFs to the corresponding master DOFs
/**
......@@ -93,7 +91,8 @@ namespace AMDiS
* columns) of dependant DOFs to unit-rows[columns] with entries of -1 at master DOF indices.
* Uses a backend-specific implementation.
**/
void apply(Mat& matrix, Sol& solution, Rhs& rhs) override;
template <class Mat, class Sol, class Rhs>
void apply(Mat& matrix, Sol& solution, Rhs& rhs);
/// Return the map of DOF associations
auto const& associations() const
......@@ -134,35 +133,35 @@ namespace AMDiS
/// Make a PeriodicBC from a subspacebasis
template <class Mat, class Sol, class Rhs, class B, class TP>
template <class B, class TP>
auto makePeriodicBC(Dune::Functions::SubspaceBasis<B, TP> basis,
BoundarySubset<typename B::GridView::Intersection> boundarySubset,
Impl::FaceTrafo<B> trafo)
{
using BcType = PeriodicBC<Mat, Sol, Rhs, B, TP>;
using BcType = PeriodicBC<B, TP>;
return BcType(std::move(basis), std::move(boundarySubset), std::move(trafo));
}
/// Make a PeriodicBC from a global basis and treepath
template <class Mat, class Sol, class Rhs, class B, class TP,
template <class B, class TP,
REQUIRES(Concepts::GlobalBasis<B>)>
auto makePeriodicBC(B const& basis, TP const& tp,
BoundarySubset<typename B::GridView::Intersection> boundarySubset,
Impl::FaceTrafo<B> trafo)
{
return makePeriodicBC<Mat, Sol, Rhs>(Dune::Functions::subspaceBasis(basis, tp),
std::move(boundarySubset), std::move(trafo));
return makePeriodicBC(Dune::Functions::subspaceBasis(basis, tp),
std::move(boundarySubset), std::move(trafo));
}
/// Make a PeriodicBC from a global basis
template <class Mat, class Sol, class Rhs, class B,
template <class B,
REQUIRES(Concepts::GlobalBasis<B>)>
auto makePeriodicBC(B const& basis,
BoundarySubset<typename B::GridView::Intersection> boundarySubset,
Impl::FaceTrafo<B> trafo)
{
return makePeriodicBC<Mat, Sol, Rhs>(basis, makeTreePath(), std::move(boundarySubset),
std::move(trafo));
return makePeriodicBC(basis, makeTreePath(), std::move(boundarySubset),
std::move(trafo));
}
} // end namespace AMDiS
......
......@@ -14,8 +14,8 @@
namespace AMDiS {
template <class Mat, class Sol, class Rhs, class Basis, class TP>
void PeriodicBC<Mat, Sol, Rhs, Basis, TP>::
template <class Basis, class TP>
void PeriodicBC<Basis, TP>::
init()
{
if (!bool(hasConnectedIntersections_)) {
......@@ -41,8 +41,8 @@ init()
}
template <class Mat, class Sol, class Rhs, class Basis, class TP>
void PeriodicBC<Mat, Sol, Rhs, Basis, TP>::
template <class Basis, class TP>
void PeriodicBC<Basis, TP>::
initAssociations()
{
using std::sqrt;
......@@ -124,8 +124,8 @@ namespace Impl
}
template <class Mat, class Sol, class Rhs, class Basis, class TP>
void PeriodicBC<Mat, Sol, Rhs, Basis, TP>::
template <class Basis, class TP>
void PeriodicBC<Basis, TP>::
initAssociations2()
{
using std::sqrt;
......@@ -191,9 +191,9 @@ initAssociations2()
}
template <class Mat, class Sol, class Rhs, class Basis, class TP>
template <class Basis, class TP>
template <class Node>
auto PeriodicBC<Mat, Sol, Rhs, Basis, TP>::
auto PeriodicBC<Basis, TP>::
coords(Node const& tree, std::vector<std::size_t> const& localIndices) const
{
std::vector<Domain> dofCoords(localIndices.size());
......@@ -231,8 +231,9 @@ coords(Node const& tree, std::vector<std::size_t> const& localIndices) const
}
template <class Mat, class Sol, class Rhs, class Basis, class TP>
void PeriodicBC<Mat, Sol, Rhs, Basis, TP>::
template <class Basis, class TP>
template <class Mat, class Sol, class Rhs>
void PeriodicBC<Basis, TP>::
apply(Mat& matrix, Sol& solution, Rhs& rhs)
{
periodicBC(matrix, solution, rhs, periodicNodes_, associations_);
......
......@@ -305,6 +305,10 @@ namespace AMDiS
void addPeriodicBC(BoundaryType id, WorldMatrix const& A, WorldVector const& b);
/** @} */
void addConstraint(BoundaryCondition<SystemMatrix, SolutionVector, SystemVector> constraint)
{
constraints_.push_back(constraint);
}
public:
......@@ -497,7 +501,8 @@ namespace AMDiS
void initGlobalBasis();
private:
protected:
/// Name of this problem.
std::string name_;
......@@ -539,8 +544,8 @@ namespace AMDiS
/// for each node in the basis tree, indexed by [to_string(treePath)][element index]
std::map<std::string, std::vector<double>> estimates_;
private: // some internal data-structures
BoundaryConditions<SystemMatrix, SolutionVector, SystemVector, GlobalBasis, GlobalBasis> boundaryConditions_;
/// List of constraints to apply to matrix, solution and rhs
std::list<BoundaryCondition<SystemMatrix, SolutionVector, SystemVector>> constraints_;
};
......
......@@ -303,17 +303,11 @@ addDirichletBC(Predicate const& predicate, RowTreePath row, ColTreePath col, Val
static_assert(isValidTreePath, "Invalid row and/or col treepath passed to addDirichletBC!");
if constexpr (isValidPredicate && isValidTreePath) {
auto localView = globalBasis_->localView();
auto i = makeTreePath(row);
auto j = makeTreePath(col);
auto rowBasis = Dune::Functions::subspaceBasis(*globalBasis_, i);
auto colBasis = Dune::Functions::subspaceBasis(*globalBasis_, j);
auto basis = Dune::Functions::subspaceBasis(*globalBasis_, makeTreePath(col));
auto valueGridFct = makeGridFunction(values, this->gridView());
auto bc = makeDirichletBC<SystemMatrix, SolutionVector, SystemVector>(
std::move(rowBasis), std::move(colBasis), {predicate}, valueGridFct);
boundaryConditions_[i][j].push_back(makeUniquePtr(std::move(bc)));
constraints_.push_back(DirichletBC{
std::move(basis), {predicate}, valueGridFct});
}
}
......@@ -330,17 +324,11 @@ addDirichletBC(BoundaryType id, RowTreePath row, ColTreePath col, Values const&
static_assert(isValidTreePath, "Invalid row and/or col treepath passed to addDirichletBC!");
if constexpr (isValidTreePath) {
auto localView = globalBasis_->localView();
auto i = makeTreePath(row);
auto j = makeTreePath(col);
auto rowBasis = Dune::Functions::subspaceBasis(*globalBasis_, i);
auto colBasis = Dune::Functions::subspaceBasis(*globalBasis_, j);
auto basis = Dune::Functions::subspaceBasis(*globalBasis_, makeTreePath(col));
auto valueGridFct = makeGridFunction(values, this->gridView());
auto bc = makeDirichletBC<SystemMatrix, SolutionVector, SystemVector>(
std::move(rowBasis), std::move(colBasis), {*boundaryManager_, id}, valueGridFct);
boundaryConditions_[i][j].push_back(makeUniquePtr(std::move(bc)));
constraints_.push_back(DirichletBC{
std::move(basis), {*boundaryManager_, id}, valueGridFct});
}
}
......@@ -351,12 +339,12 @@ addPeriodicBC(BoundaryType id, WorldMatrix const& matrix, WorldVector const& vec
{
auto localView = globalBasis_->localView();
auto basis = Dune::Functions::subspaceBasis(*globalBasis_, makeTreePath());
auto bc = makePeriodicBC<SystemMatrix, SolutionVector, SystemVector>(
std::move(basis), {*boundaryManager_, id}, {matrix, vector});
boundaryConditions_[makeTreePath()][makeTreePath()].push_back(makeUniquePtr(std::move(bc)));
constraints_.push_back(makePeriodicBC(
std::move(basis), {*boundaryManager_, id}, {matrix, vector}));
}
template <class Traits>
void ProblemStat<Traits>::
solve(AdaptInfo& /*adaptInfo*/, bool createMatrixData, bool storeMatrixData)
......@@ -470,13 +458,9 @@ buildAfterAdapt(AdaptInfo& /*adaptInfo*/, Flag /*flag*/, bool asmMatrix, bool as
Dune::Timer t;
Dune::Timer t2;
auto localView = globalBasis_->localView();
for_each_node(localView.tree(), [&](auto&&, auto rowTp) -> void {
for_each_node(localView.tree(), [&](auto&&, auto colTp) -> void {
for (auto bc : boundaryConditions_[rowTp][colTp])
bc->init();
});
});
// 0. initialize boundary condition and other constraints
for (auto& bc : constraints_)
bc.init();
t2.reset();
......@@ -491,6 +475,7 @@ buildAfterAdapt(AdaptInfo& /*adaptInfo*/, Flag /*flag*/, bool asmMatrix, bool as
msg("{} local DOFs", rhs_->localSize());
// 2. traverse grid and assemble operators on the elements
auto localView = globalBasis_->localView();
for (auto const& element : elements(gridView(), PartitionSet{})) {
localView.bind(element);
......@@ -510,13 +495,10 @@ buildAfterAdapt(AdaptInfo& /*adaptInfo*/, Flag /*flag*/, bool asmMatrix, bool as
t2.reset();
solution_->resize(sizeInfo(*globalBasis_));
for_each_node(localView.tree(), [&](auto&&, auto rowTp) -> void {