version 0.4.1
Ikarus::DirichletValues< B, FC > Class Template Reference

Class for handling Dirichlet boundary conditions in Ikarus. More...

#include <ikarus/utils/dirichletvalues.hh>

Public Types

using Basis = std::remove_cvref_t< B >
 
using FlagsType = FC
 
using BackendType = decltype(Dune::Functions::istlVectorBackend(std::declval< FlagsType & >()))
 

Public Member Functions

 DirichletValues (const B &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
 

Detailed Description

template<typename B, typename FC = std::vector<bool>>
class Ikarus::DirichletValues< B, FC >

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.

Template Parameters
BType of the finite element basis
FCType for storing Dirichlet flags (default is std::vector<bool>)

Member Typedef Documentation

◆ BackendType

template<typename B , typename FC = std::vector<bool>>
using Ikarus::DirichletValues< B, FC >::BackendType = decltype(Dune::Functions::istlVectorBackend(std::declval<FlagsType&>()))

◆ Basis

template<typename B , typename FC = std::vector<bool>>
using Ikarus::DirichletValues< B, FC >::Basis = std::remove_cvref_t<B>

◆ FlagsType

template<typename B , typename FC = std::vector<bool>>
using Ikarus::DirichletValues< B, FC >::FlagsType = FC

Constructor & Destructor Documentation

◆ DirichletValues()

template<typename B , typename FC = std::vector<bool>>
Ikarus::DirichletValues< B, FC >::DirichletValues ( const B &  basis)
inlineexplicit

Member Function Documentation

◆ basis()

template<typename B , typename FC = std::vector<bool>>
const auto & Ikarus::DirichletValues< B, FC >::basis ( ) const
inline

◆ container()

template<typename B , typename FC = std::vector<bool>>
auto & Ikarus::DirichletValues< B, FC >::container ( ) const
inline

◆ evaluateInhomogeneousBoundaryCondition()

template<typename B , typename FC = std::vector<bool>>
void Ikarus::DirichletValues< B, FC >::evaluateInhomogeneousBoundaryCondition ( Eigen::VectorXd &  xIh,
const double &  lambda 
)
inline

This function evaluates all stored inhomogeneous Dirichlet boundary functions at all positions where the corresponding degrees of freedom are true.

Parameters
xIhThe vector where the interpolated result should be stored
lambdaThe load factor

◆ evaluateInhomogeneousBoundaryConditionDerivative()

template<typename B , typename FC = std::vector<bool>>
void Ikarus::DirichletValues< B, FC >::evaluateInhomogeneousBoundaryConditionDerivative ( Eigen::VectorXd &  xIh,
const double &  lambda 
)
inline

This function evaluates all stored inhomogeneous Dirichlet boundary derivative functions at all positions where the corresponding degrees of freedom are true.

Parameters
xIhThe vector where the interpolated result should be stored
lambdaThe load factor

◆ fixBoundaryDOFs()

template<typename B , typename FC = std::vector<bool>>
template<typename F >
void Ikarus::DirichletValues< B, FC >::fixBoundaryDOFs ( F &&  f)
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.

Parameters
fA callback function
Here is the caller graph for this function:

◆ fixDOFs()

template<typename B , typename FC = std::vector<bool>>
template<typename F >
void Ikarus::DirichletValues< B, FC >::fixDOFs ( F &&  f)
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.

Parameters
fA callback function
Here is the caller graph for this function:

◆ fixedDOFsize()

template<typename B , typename FC = std::vector<bool>>
auto Ikarus::DirichletValues< B, FC >::fixedDOFsize ( ) const
inline

◆ fixIthDOF()

template<typename B , typename FC = std::vector<bool>>
void Ikarus::DirichletValues< B, FC >::fixIthDOF ( typename Basis::MultiIndex  i)
inline
Parameters
iAn index indicating the DOF number to be fixed

◆ isConstrained() [1/2]

template<typename B , typename FC = std::vector<bool>>
template<typename MultiIndex >
requires (not std::integral<MultiIndex>)
bool Ikarus::DirichletValues< B, FC >::isConstrained ( const MultiIndex &  multiIndex) const
inline

◆ isConstrained() [2/2]

template<typename B , typename FC = std::vector<bool>>
bool Ikarus::DirichletValues< B, FC >::isConstrained ( std::size_t  i) const
inline

◆ size()

template<typename B , typename FC = std::vector<bool>>
auto Ikarus::DirichletValues< B, FC >::size ( ) const
inline

◆ storeInhomogeneousBoundaryCondition()

template<typename B , typename FC = std::vector<bool>>
template<typename F >
void Ikarus::DirichletValues< B, FC >::storeInhomogeneousBoundaryCondition ( F &&  f)
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.

Parameters
fA callback function

Member Data Documentation

◆ worldDimension

template<typename B , typename FC = std::vector<bool>>
constexpr int Ikarus::DirichletValues< B, FC >::worldDimension = Basis::GridView::dimensionworld
staticconstexpr

The documentation for this class was generated from the following file: