19#include <dune/common/float_cmp.hh>
20#include <dune/common/hybridutilities.hh>
21#include <dune/common/indices.hh>
22#include <dune/functions/backends/istlvectorbackend.hh>
23#include <dune/functions/functionspacebases/boundarydofs.hh>
24#include <dune/functions/functionspacebases/compositebasis.hh>
25#include <dune/functions/functionspacebases/flatmultiindex.hh>
26#include <dune/functions/functionspacebases/interpolate.hh>
27#include <dune/functions/functionspacebases/subspacebasis.hh>
31#include <autodiff/forward/real/real.hpp>
37template <
class K,
int N>
49template <Concepts::EigenVector T>
73template <
typename B,
typename FC = std::vector<
bool>>
77 using Basis = std::remove_cvref_t<B>;
79 using BackendType =
decltype(Dune::Functions::istlVectorBackend(std::declval<FlagsType&>()));
83 using Tree = LocalView::Tree;
88 dirichletFlagsBackend_{dirichletFlags_} {
89 dirichletFlagsBackend_.resize(basis_);
90 std::fill(dirichletFlags_.begin(), dirichletFlags_.end(),
false);
91 if constexpr (Tree::isComposite) {
92 Dune::Hybrid::forEach(
93 Dune::Hybrid::integralRange(Dune::index_constant<std::tuple_size_v<typename Basis::PreBasis::SubPreBases>>()),
95 static_assert(not Tree::template Child<i>::Type::isComposite,
96 "DirichletValues is not implemented to handle a composite basis within a composite basis.");
110 template <
typename F,
typename TreePath = Dune::TypeTree::Hybr
idTreePath<>>
112 using namespace Dune::Functions;
113 using SubSpaceLocalView =
114 typename std::remove_cvref_t<decltype(subspaceBasis(basis_, std::forward<TreePath>(treePath)))>
::LocalView;
117 auto lambda = [&](
auto&& indexGlobal) { f(dirichletFlagsBackend_, indexGlobal); };
118 Dune::Functions::forEachBoundaryDOF(subspaceBasis(basis_, std::forward<TreePath>(treePath)), lambda);
119 }
else if constexpr (Concepts::IsFunctorWithArgs<F, BackendType, int, SubSpaceLocalView>) {
120 auto lambda = [&](
auto&& localIndex,
auto&& localView) { f(dirichletFlagsBackend_, localIndex, localView); };
121 Dune::Functions::forEachBoundaryDOF(subspaceBasis(basis_, std::forward<TreePath>(treePath)), lambda);
122 }
else if constexpr (Concepts::IsFunctorWithArgs<F,
BackendType, int, SubSpaceLocalView,
123 typename Basis::GridView::Intersection>) {
124 auto lambda = [&](
auto&& localIndex,
auto&& localView,
auto&& intersection) {
125 f(dirichletFlagsBackend_, localIndex, localView, intersection);
127 Dune::Functions::forEachBoundaryDOF(subspaceBasis(basis_, std::forward<TreePath>(treePath)), lambda);
129 static_assert(Dune::AlwaysFalse<F>(),
"fixBoundaryDOFs: A function with this signature is not supported");
141 template <
typename F>
143 f(basis_, dirichletFlagsBackend_);
152 template <
typename MultiIndex>
153 requires(not std::integral<MultiIndex>)
155 dirichletFlagsBackend_[i] = flag;
166 requires(std::same_as<typename Basis::MultiIndex, Dune::Functions::FlatMultiIndex<size_t>>)
168 dirichletFlags_[i] = flag;
174 void reset() { std::fill(dirichletFlags_.begin(), dirichletFlags_.end(),
false); }
177 const auto&
basis()
const {
return basis_; }
180 template <
typename MultiIndex>
181 requires(not std::integral<MultiIndex>)
183 return dirichletFlagsBackend_[multiIndex];
188 requires(std::same_as<typename Basis::MultiIndex, Dune::Functions::FlatMultiIndex<size_t>>)
190 return dirichletFlags_[i];
194 auto fixedDOFsize()
const {
return std::ranges::count(dirichletFlags_,
true); }
197 auto size()
const {
return dirichletFlags_.size(); }
213 template <
typename F>
215 auto derivativeLambda = [&](
const auto& globalCoord,
const double& lambda) {
216 autodiff::real lambdaDual = lambda;
218 return derivative(f(globalCoord, lambdaDual));
220 dirichletFunctions_.push_back({f, derivativeLambda});
221 setInhomogeneousBoundaryConditionFlag(lambda);
230 for (Eigen::Index i = 0; i < xIh.size(); ++i)
245 Eigen::VectorXd inhomogeneousBoundaryVectorDummy;
246 inhomogeneousBoundaryVectorDummy.setZero(this->
size());
247 xIh.resizeLike(inhomogeneousBoundaryVectorDummy);
249 auto runInterpolateImpl = [&](
const auto&
basis) {
250 for (
auto& f : dirichletFunctions_) {
251 interpolate(
basis, inhomogeneousBoundaryVectorDummy,
252 [&](
const auto& globalCoord) {
return f.value(globalCoord, lambda); });
253 xIh += inhomogeneousBoundaryVectorDummy;
256 runInterpolate(runInterpolateImpl);
269 Eigen::VectorXd inhomogeneousBoundaryVectorDummy;
270 inhomogeneousBoundaryVectorDummy.setZero(this->
size());
271 xIh.resizeLike(inhomogeneousBoundaryVectorDummy);
273 auto runInterpolateImpl = [&](
const auto&
basis) {
274 for (
auto& f : dirichletFunctions_) {
275 interpolate(
basis, inhomogeneousBoundaryVectorDummy,
276 [&](
const auto& globalCoord) {
return f.derivative(globalCoord, lambda); });
277 xIh += inhomogeneousBoundaryVectorDummy;
280 runInterpolate(runInterpolateImpl);
287 struct DirichletFunctions
291 Signature derivative;
293 std::vector<DirichletFunctions> dirichletFunctions_;
295 void setInhomogeneousBoundaryConditionFlag(
double lambda) {
296 Eigen::VectorXd inhomogeneousBoundaryVectorDummy;
297 inhomogeneousBoundaryVectorDummy.setZero(this->
size());
299 for (
const std::size_t i : Dune::range(this->
size()))
300 if (Dune::FloatCmp::ne(inhomogeneousBoundaryVectorDummy[i], 0.0))
304 template <
typename F>
305 void runInterpolate(F&& f)
const {
306 if constexpr (Tree::isComposite)
307 f(subspaceBasis(basis_, Dune::Indices::_0));
Implementation of creating a flat basis from a possibly blocked basis.
Definition: assemblermanipulatorbuildingblocks.hh:22
Definition: utils/dirichletvalues.hh:36
Definition: utils/dirichletvalues.hh:38
A helper struct to derive the SizeType of the underlying container.
Definition: utils/dirichletvalues.hh:47
Eigen::Index SizeType
Definition: utils/dirichletvalues.hh:52
std::vector< bool >::size_type SizeType
Definition: utils/dirichletvalues.hh:58
Class for handling Dirichlet boundary conditions in Ikarus.
Definition: utils/dirichletvalues.hh:75
auto size() const
Definition: utils/dirichletvalues.hh:197
void setZeroAtConstrainedDofs(Eigen::VectorXd &xIh) const
Function to write zeros at constrained Dirichlet entries.
Definition: utils/dirichletvalues.hh:229
void setSingleDOF(const MultiIndex i, bool flag)
Fixes and unfixes (set boolean value to true or false) a specific degree of freedom.
Definition: utils/dirichletvalues.hh:154
bool isConstrained(const MultiIndex &multiIndex) const
Definition: utils/dirichletvalues.hh:182
bool isConstrained(std::size_t i) const
Definition: utils/dirichletvalues.hh:187
void setSingleDOF(std::size_t i, bool flag)
Fixes or unfixes (set boolean value to true or false) a specific degree of freedom.
Definition: utils/dirichletvalues.hh:165
typename utils::Impl::PreBasisInfo< Tree >::NodalSolutionType NodalSolutionType
Definition: utils/dirichletvalues.hh:84
std::remove_cvref_t< B > Basis
Definition: utils/dirichletvalues.hh:77
void fixBoundaryDOFs(F &&f, TreePath &&treePath={})
Function to fix (set boolean values to true or false) degrees of freedom on the boundary.
Definition: utils/dirichletvalues.hh:111
decltype(Dune::Functions::istlVectorBackend(std::declval< FlagsType & >())) BackendType
Definition: utils/dirichletvalues.hh:79
void evaluateInhomogeneousBoundaryConditionDerivative(Eigen::VectorXd &xIh, const double &lambda) const
Function to evaluate all stored inhomogeneous Dirichlet boundary derivative functions.
Definition: utils/dirichletvalues.hh:268
void storeInhomogeneousBoundaryCondition(F &&f, double lambda=1.0)
Function to insert a function for inhomogeneous Dirichlet boundary conditions.
Definition: utils/dirichletvalues.hh:214
LocalView::Tree Tree
Definition: utils/dirichletvalues.hh:83
static constexpr int worldDimension
Definition: utils/dirichletvalues.hh:81
void fixDOFs(F &&f)
Function to fix (set boolean values to true or false) degrees of freedom.
Definition: utils/dirichletvalues.hh:142
FC FlagsType
Definition: utils/dirichletvalues.hh:78
void evaluateInhomogeneousBoundaryCondition(Eigen::VectorXd &xIh, const double &lambda) const
Function to evaluate all stored inhomogeneous Dirichlet boundary functions.
Definition: utils/dirichletvalues.hh:244
Basis::LocalView LocalView
Definition: utils/dirichletvalues.hh:82
void reset()
Resets all degrees of freedom.
Definition: utils/dirichletvalues.hh:174
typename DeriveSizeType< FlagsType >::SizeType SizeType
Definition: utils/dirichletvalues.hh:80
const auto & basis() const
Definition: utils/dirichletvalues.hh:177
auto fixedDOFsize() const
Definition: utils/dirichletvalues.hh:194
DirichletValues(const B &basis)
Definition: utils/dirichletvalues.hh:86
static constexpr std::size_t numberOfChildrenAtNode
Definition: utils/dirichletvalues.hh:85
auto & container() const
Definition: utils/dirichletvalues.hh:200
Concept defining the requirements for functors with arguments.
Definition: utils/concepts.hh:357