7 #include <amdis/Boundary.hpp> 8 #include <amdis/BoundaryCondition.hpp> 9 #include <amdis/BoundarySubset.hpp> 10 #include <amdis/common/Concepts.hpp> 11 #include <amdis/common/TypeTraits.hpp> 12 #include <amdis/typetree/RangeType.hpp> 13 #include <amdis/typetree/TreePath.hpp> 36 template <
class Basis,
class RowPath,
class ColPath,
class ValueGr
idFct>
39 using GridView =
typename Basis::GridView;
40 using Intersection =
typename GridView::Intersection;
42 using Domain =
typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
43 using Range = TYPEOF(std::declval<ValueGridFct>()(std::declval<Domain>()));
48 DirichletBC(B&& basis, RowPath
const& row, ColPath
const& col,
50 : basis_{wrap_or_share(FWD(basis))}
53 , boundarySubset_(std::move(boundarySubset))
54 , valueGridFct_(std::move(values))
60 :
DirichletBC(FWD(basis), makeTreePath(), makeTreePath(),
61 std::move(boundarySubset),
std::move(values))
78 template <
class Mat,
class Sol,
class Rhs>
79 void apply(Mat& matrix, Sol& solution, Rhs& rhs,
bool symmetric =
false);
82 std::shared_ptr<Basis const> basis_;
86 ValueGridFct valueGridFct_;
88 std::vector<typename Basis::MultiIndex> rowIndices_;
89 std::vector<typename Basis::MultiIndex> colIndices_;
95 template <
class B,
class Row,
class Col,
class Values,
class Basis = Underlying_t<B>,
96 REQUIRES(Concepts::GlobalBasis<Basis>)>
101 template <
class B,
class Values,
class Basis = Underlying_t<B>,
102 REQUIRES(Concepts::GlobalBasis<Basis>)>
108 #include "DirichletBC.inc.hpp"
DirichletBC(B &&basis, BoundarySubset< Intersection > boundarySubset, ValueGridFct values)
Make a DirichletBC from a global basis.
Definition: DirichletBC.hpp:59
Definition: FieldMatVec.hpp:12
Implements a boundary condition of Dirichlet-type.
Definition: DirichletBC.hpp:37
Definition: AdaptBase.hpp:6
DirichletBC(B &&basis, RowPath const &row, ColPath const &col, BoundarySubset< Intersection > boundarySubset, ValueGridFct values)
Make a DirichletBC from a basis with treepath arguments.
Definition: DirichletBC.hpp:48
void init()
Definition: DirichletBC.inc.hpp:18
void apply(Mat &matrix, Sol &solution, Rhs &rhs, bool symmetric=false)
Apply dirichlet BC to matrix and vector.
Definition: DirichletBC.inc.hpp:81