Class for handling Dirichlet boundary conditions in Ikarus. More...
#include <ikarus/utils/dirichletvalues.hh>
Public Types | |
using | Basis = std::remove_cvref_t< Basis_ > |
using | FlagsType = FlagsType_ |
using | BackendType = decltype(Dune::Functions::istlVectorBackend(std::declval< FlagsType & >())) |
Public Member Functions | |
DirichletValues (const Basis_ &p_basis) | |
template<typename F > | |
void | fixBoundaryDOFs (F &&f) |
Function to fix (set boolean values to true or false) degrees of freedom on the boundary. More... | |
template<typename F > | |
void | fixDOFs (F &&f) |
Function to fix (set boolean values to true or false) degrees of freedom. More... | |
void | fixIthDOF (typename Basis::MultiIndex i) |
Function to fix (set boolean values to true or false) of degrees of freedom. More... | |
const auto & | basis () const |
template<typename MultiIndex > requires (not std::integral<MultiIndex>) | |
bool | isConstrained (const MultiIndex &multiIndex) const |
bool | isConstrained (std::size_t i) const |
auto | fixedDOFsize () const |
auto | size () const |
auto & | container () const |
template<typename F > | |
void | storeInhomogeneousBoundaryCondition (F &&f) |
Function to insert a function for inhomogeneous Dirichlet boundary conditions. More... | |
void | evaluateInhomogeneousBoundaryCondition (Eigen::VectorXd &xIh, const double &lambda) |
Function to evaluate all stored inhomogeneous Dirichlet boundary functions. More... | |
void | evaluateInhomogeneousBoundaryConditionDerivative (Eigen::VectorXd &xIh, const double &lambda) |
Function to evaluate all stored inhomogeneous Dirichlet boundary derivative functions. More... | |
Static Public Attributes | |
static constexpr int | worldDimension = Basis::GridView::dimensionworld |
The DirichletValues class provides functionalities for fixing degrees of freedom and storing inhomogeneous Dirichlet boundary conditions. It supports fixing degrees of freedom using various callback functions and stores functions for inhomogeneous Dirichlet boundary conditions.
Basis_ | Type of the finite element basis |
FlagsType_ | Type for storing Dirichlet flags (default is std::vector<bool>) |
using Ikarus::DirichletValues< Basis_, FlagsType_ >::BackendType = decltype(Dune::Functions::istlVectorBackend(std::declval<FlagsType&>())) |
using Ikarus::DirichletValues< Basis_, FlagsType_ >::Basis = std::remove_cvref_t<Basis_> |
using Ikarus::DirichletValues< Basis_, FlagsType_ >::FlagsType = FlagsType_ |
|
inlineexplicit |
|
inline |
|
inline |
|
inline |
This function evaluates all stored inhomogeneous Dirichlet boundary functions at all positions where the corresponding degrees of freedom are true.
xIh | The vector where the interpolated result should be stored |
lambda | The load factor |
|
inline |
This function evaluates all stored inhomogeneous Dirichlet boundary derivative functions at all positions where the corresponding degrees of freedom are true.
xIh | The vector where the interpolated result should be stored |
lambda | The load factor |
|
inline |
This function takes a callback function, which will be called with the boolean vector of fixed boundary degrees of freedom and the usual arguments of Dune::Functions::forEachBoundaryDOF
.
f | A callback function |
|
inline |
This function takes a callback function, which will be called with the stored function space basis and the boolean vector of fixed boundary degrees of freedom.
f | A callback function |
|
inline |
|
inline |
i | An index indicating the DOF number to be fixed |
|
inline |
|
inline |
|
inline |
|
inline |
This function takes a callback function, which will be called with the current coordinate vector and the scalar load factor. It creates internally the first derivative of the passed function and stores them simultaneously.
f | A callback function |
|
staticconstexpr |