Skip to content
Snippets Groups Projects

Allow storage of GridFunctions in DirichletBC

Merged Praetorius, Simon requested to merge feature/dirichletbc_with_gridfunction into master
2 files
+ 25
32
Compare changes
  • Side-by-side
  • Inline
Files
2
+ 19
24
@@ -40,17 +40,18 @@ namespace AMDiS
* \tparam CB Basis of the column FE space
* \tparam CTP Treepath to the column node this constraint applies to
**/
template <class Mat, class Sol, class Rhs, class RB, class RTP, class CB, class CTP>
template <class Mat, class Sol, class Rhs, class RowSubBasis, class ColSubBasis, class ValueGridFct>
class DirichletBC
: public BoundaryCondition<Mat, Sol, Rhs>
{
using Super = BoundaryCondition<Mat, Sol, Rhs>;
using RowSubBasis = Dune::Functions::SubspaceBasis<RB, RTP>;
using ColSubBasis = Dune::Functions::SubspaceBasis<CB, CTP>;
// We assume CB's GridView to be the same as RB's
using Domain = typename RB::GridView::template Codim<0>::Geometry::GlobalCoordinate;
using Intersection = typename RB::GridView::Intersection;
using Range = RangeType_t<typename ColSubBasis::LocalView::Tree>;
using GridView = typename RowSubBasis::GridView;
using Intersection = typename GridView::Intersection;
using Domain = typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
using Range = TYPEOF(std::declval<ValueGridFct>()(std::declval<Domain>()));
public:
/// Constructor accepting subspacebases.
@@ -64,17 +65,6 @@ namespace AMDiS
, values_(FWD(values))
{}
/// Constructor accepting global bases and treepaths.
template <class Values,
REQUIRES(Concepts::Functor<Values, Range(Domain)>)>
DirichletBC(RB const& rowBasis, RTP const& rowTreePath,
CB const& colBasis, CTP const& colTreePath,
BoundarySubset<Intersection> boundarySubset, Values&& values)
: DirichletBC(Dune::Functions::subspaceBasis(rowBasis, rowTreePath),
Dune::Functions::subspaceBasis(colBasis, colTreePath),
std::move(boundarySubset), FWD(values))
{}
// fill \ref dirichletNodes_ with 1 or 0 whether DOF corresponds to the boundary or not.
/**
* \see BoundaryCondition::init
@@ -92,7 +82,8 @@ namespace AMDiS
RowSubBasis rowBasis_;
ColSubBasis colBasis_;
BoundarySubset<Intersection> boundarySubset_;
std::function<Range(Domain)> values_;
ValueGridFct values_;
std::vector<bool> dirichletNodes_;
};
@@ -104,7 +95,9 @@ namespace AMDiS
BoundarySubset<typename RB::GridView::Intersection> boundarySubset,
Values&& values)
{
using BcType = DirichletBC<Mat, Sol, Rhs, RB, RTP, CB, CTP>;
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));
}
@@ -117,9 +110,10 @@ namespace AMDiS
BoundarySubset<typename RB::GridView::Intersection> boundarySubset,
Values&& values)
{
using BcType = DirichletBC<Mat, Sol, Rhs, RB, RTP, CB, CTP>;
return BcType(rowBasis, rowTreePath, colBasis, colTreePath,
std::move(boundarySubset), FWD(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
@@ -130,8 +124,9 @@ namespace AMDiS
BoundarySubset<typename RB::GridView::Intersection> boundarySubset,
Values&& values)
{
return makeDirichletBC<Mat, Sol, Rhs>(rowBasis, makeTreePath(), colBasis, makeTreePath(),
std::move(boundarySubset), FWD(values));
return makeDirichletBC<Mat, Sol, Rhs>(
rowBasis, makeTreePath(), colBasis, makeTreePath(),
std::move(boundarySubset), FWD(values));
}
} // end namespace AMDiS
Loading