version 0.4.1

Collection of several utilities. More...

Collaboration diagram for Utilities:

Classes

class  Ikarus::BasisHandler< PB >
 Wrapper class for a hierarchical basis constructed from a pre-basis. More...
 
class  Ikarus::DirichletValues< B, FC >
 Class for handling Dirichlet boundary conditions in Ikarus. More...
 
struct  Ikarus::FlatPreBasis< PreBasis >
 Transform a PreBasis into one with flat index-merging strategyThis utility takes a pre-basis and converts recursively all index-merging strategies into their flat analog, i.e. BlockedInterleaved is converted into FlatInterleaved and BlockedLexicographic is transformed into FlatLexicographic. More...
 
class  Ikarus::NonLinearOperator< TypeListOne, TypeListTwo >
 Represents a NonLinearOperator class for handling nonlinear operators. More...
 
struct  Ikarus::utils::SolverDefault
 Default functor for solving operations. More...
 
struct  Ikarus::utils::UpdateDefault
 Default functor for updating operations. More...
 

Macros

#define MAKE_ENUM(type, ...)
 Macro to create an enumeration with a string conversion function.The macro creates an enum class with a BEGIN and END enumerator, and provides a constexpr toString function for string conversion. More...
 

Functions

template<typename Fun , typename... Vars, typename... Args, typename U , typename G , typename H >
void Ikarus::utils::hessianN (const Fun &f, const autodiff::Wrt< Vars... > &wrt, const autodiff::At< Args... > &at, U &u, std::array< G, U::RowsAtCompileTime > &g, std::array< H, U::RowsAtCompileTime > &h)
 Computes the Hessian matrix for each parameter of a given function.The Hessian matrix represents the second-order partial derivatives of the function with respect to the specified variables. More...
 
std::tuple< Dune::Functions::Polynomial< double >, decltype(Eigen::seq(0, 0))> Ikarus::utils::findLineSegment (const Eigen::VectorXd &x, const Eigen::VectorXd &y, int segmentSize)
 Find a linear segment in a set of data points. More...
 
template<class PreBasis >
decltype(auto) Ikarus::flatPreBasis (const PreBasis &preBasis)
 Generator function for a flatted PreBasis. More...
 
template<int size, typename LV >
void Ikarus::utils::obtainLagrangeNodePositions (const LV &localView, std::vector< Dune::FieldVector< double, size > > &lagrangeNodeCoords)
 A function to obtain the global positions of the nodes of an element with Lagrangian basis, see Dune book page 314. More...
 
template<typename FE >
auto Ikarus::utils::referenceElementSubEntityPositions (FE &fe, int codim)
 A function to obtain the local coordinates of subentities of an FiniteElement. More...
 
template<typename FE >
auto Ikarus::utils::referenceElementVertexPositions (FE &fe)
 A function to obtain the local coordinates the vertices of an FiniteElement. More...
 
template<typename NonlinearOperator , typename UpdateType = typename NonlinearOperator::template ParameterValue<0>>
bool Ikarus::utils::checkGradient (NonlinearOperator &nonLinOp, CheckFlags checkFlags=CheckFlags(), std::function< void(typename NonlinearOperator::template ParameterValue< 0 > &, const UpdateType &)> p_updateFunction=[](typename NonlinearOperator::template ParameterValue< 0 > &a, const UpdateType &b) { a+=b;})
 Checks the gradient of a nonlinear operator. More...
 
template<typename NonlinearOperator , typename UpdateType = typename NonlinearOperator::template ParameterValue<0>>
bool Ikarus::utils::checkJacobian (NonlinearOperator &nonLinOp, CheckFlags checkFlags=CheckFlags(), std::function< void(typename NonlinearOperator::template ParameterValue< 0 > &, const UpdateType &)> p_updateFunction=[](typename NonlinearOperator::template ParameterValue< 0 > &a, const UpdateType &b) { a+=b;})
 Checks the Jacobian of a nonlinear operator. More...
 
template<typename NonlinearOperator , typename UpdateType = typename NonlinearOperator::template ParameterValue<0>>
bool Ikarus::utils::checkHessian (NonlinearOperator &nonLinOp, CheckFlags checkFlags=CheckFlags(), std::function< void(typename NonlinearOperator::template ParameterValue< 0 > &, const UpdateType &)> p_updateFunction=[](typename NonlinearOperator::template ParameterValue< 0 > &a, const UpdateType &b) { a+=b;})
 Checks the Hessian of a nonlinear operator. More...
 
template<typename Derived >
auto Ikarus::orthonormalizeMatrixColumns (const Eigen::MatrixBase< Derived > &A)
 Orthonormalizes all Matrix columns using Gram-Schmidt Orthogonalization. More...
 
template<typename ValueType >
auto Ikarus::viewAsFlatEigenVector (Dune::BlockVector< ValueType > &blockedVector)
 View Dune::BlockVector as an Eigen::Vector. More...
 
template<typename ValueType >
auto Ikarus::viewAsFlatEigenVector (const Dune::BlockVector< ValueType > &blockedVector)
 View const Dune::BlockVector as an Eigen::Vector. More...
 
template<typename ValueType >
auto Ikarus::viewAsEigenMatrixAsDynFixed (Dune::BlockVector< ValueType > &blockedVector)
 View Dune::BlockVector as an Eigen::Matrix with dynamic rows and fixed columns depending on the size of the ValueType. More...
 
template<typename ValueType >
auto Ikarus::viewAsEigenMatrixAsDynFixed (const Dune::BlockVector< ValueType > &blockedVector)
 Const view Dune::BlockVector as an Eigen::Matrix with dynamic rows and fixed columns depending on the size of the ValueType. More...
 
template<typename ValueType >
auto Ikarus::viewAsEigenMatrixFixedDyn (Dune::BlockVector< ValueType > &blockedVector)
 View Dune::BlockVector as an Eigen::Matrix with fixed rows depending on the size of the ValueType and dynamic columns. More...
 
template<typename ValueType >
auto Ikarus::viewAsEigenMatrixFixedDyn (const Dune::BlockVector< ValueType > &blockedVector)
 Const view Dune::BlockVector as an Eigen::Matrix with fixed rows depending on the size of the ValueType and dynamic columns. More...
 
template<typename Type >
requires requires { Type::correctionSize; }
size_t Ikarus::correctionSize (const Dune::BlockVector< Type > &a)
 Returns the total correction size of a block vector with a Manifold as the underlying type. More...
 
template<typename Type >
requires requires { Type::valueSize; }
size_t Ikarus::valueSize (const Dune::BlockVector< Type > &a)
 Returns the total value size of a block vector with a Manifold as the underlying type. More...
 
template<typename Type , typename Derived >
requires (Ikarus::Concepts::AddAssignAble<Type, decltype(b.template segment<Type::correctionSize>(0))> and requires() { Type::correctionSize; })
Dune::BlockVector< Type > & Ikarus::operator+= (Dune::BlockVector< Type > &a, const Eigen::MatrixBase< Derived > &b)
 Enables the += operator for Dune::BlockVector += Eigen::Vector. More...
 
template<typename Type , typename Derived >
requires (Ikarus::Concepts::AddAssignAble<Type, decltype(b.template segment<Type::correctionSize>(0))> and requires() { Type::correctionSize; })
Dune::BlockVector< Type > & Ikarus::operator-= (Dune::BlockVector< Type > &a, const Eigen::MatrixBase< Derived > &b)
 Enables the -= operator for Dune::BlockVector += Eigen::Vector. More...
 
template<typename... Types, typename Derived >
Dune::TupleVector< Types... > & Ikarus::operator+= (Dune::TupleVector< Types... > &a, const Eigen::MatrixBase< Derived > &b)
 Enables the += operator for Dune::TupleVector += Eigen::Vector. More...
 
template<typename ManifoldPoint , typename Derived >
requires (Ikarus::Concepts::AddAssignAble<ManifoldPoint, decltype(b.template segment<ManifoldPoint::valueSize>(0))> and requires() { ManifoldPoint::valueSize; })
Dune::BlockVector< ManifoldPoint > & Ikarus::addInEmbedding (Dune::BlockVector< ManifoldPoint > &a, const Eigen::MatrixBase< Derived > &b)
 Enables the addition in the embedding space of a vector in the space M^n, where M is a manifold with the points of type ManifoldPoint. More...
 
template<typename Derived >
requires (!std::floating_point<Derived>)
auto Ikarus::norm (const Eigen::MatrixBase< Derived > &v)
 Adding free norm function to Eigen types. More...
 
auto Ikarus::norm (const std::floating_point auto &v)
 Helper Free Function to have the same interface as for Eigen Vector Types. More...
 
template<typename Scalar , int size>
auto Ikarus::operator* (const Eigen::DiagonalMatrix< Scalar, size > &a, const Eigen::DiagonalMatrix< Scalar, size > &b)
 Eigen::DiagonalMatrix Product Missing in Eigen. More...
 
template<typename Scalar , int size>
auto Ikarus::operator+= (Eigen::DiagonalMatrix< Scalar, size > &a, const Eigen::DiagonalMatrix< Scalar, size > &b)
 In-place addition for Eigen::DiagonalMatrix. More...
 
template<typename Derived , typename Scalar , int size>
auto Ikarus::operator+ (const Eigen::MatrixBase< Derived > &a, const Eigen::DiagonalMatrix< Scalar, size > &b)
 Eigen::Matrix + Eigen::DiagonalMatrix addition missing in Eigen. More...
 
template<typename Derived , typename Scalar , int size>
auto Ikarus::operator+ (const Eigen::DiagonalMatrix< Scalar, size > &a, const Eigen::MatrixBase< Derived > &b)
 Eigen::DiagonalMatrix + Eigen::Matrix addition missing in Eigen. More...
 
template<typename Scalar , int size>
auto Ikarus::operator- (const Eigen::DiagonalMatrix< Scalar, size > &a)
 Unary minus for Eigen::DiagonalMatrix. More...
 
template<typename Derived , typename Derived2 >
auto Ikarus::operator+ (const Eigen::MatrixBase< Derived > &a, const Eigen::DiagonalWrapper< Derived2 > &b)
 Addition of Eigen::Matrix and Eigen::DiagonalWrapper. More...
 
template<typename Derived , typename Derived2 >
auto Ikarus::operator+ (const Eigen::DiagonalWrapper< Derived > &a, const Eigen::MatrixBase< Derived2 > &b)
 Addition of Eigen::DiagonalWrapper and Eigen::Matrix. More...
 
template<typename Scalar , int size>
std::ostream & Ikarus::operator<< (std::ostream &os, const Eigen::DiagonalMatrix< Scalar, size > &a)
 Output stream operator for Eigen::DiagonalMatrix. More...
 
template<typename Derived >
Derived Ikarus::sym (const Eigen::MatrixBase< Derived > &A)
 Returns the symmetric part of a matrix. More...
 
template<typename Derived >
Derived Ikarus::skew (const Eigen::MatrixBase< Derived > &A)
 Returns the skew part of a matrix. More...
 
template<typename Derived >
void Ikarus::printForMaple (const Eigen::EigenBase< Derived > &A)
 Method to print the matrix in a format that can directly be copied to Maple. More...
 
template<typename FieldVectorT >
auto Ikarus::createRandomVector (typename FieldVectorT::value_type lower=-1, typename FieldVectorT::value_type upper=1)
 Creates a random vector of the specified type within a given range. More...
 
template<typename ScalarType >
Eigen::Matrix< ScalarType, 3, 3 > Ikarus::skew (const Eigen::Vector< ScalarType, 3 > &a)
 Create skew 3x3 matrix from 3d vector. More...
 
template<typename Derived , size_t sizeOfCondensedIndices>
auto Ikarus::staticCondensation (const Eigen::MatrixBase< Derived > &E, const std::array< size_t, sizeOfCondensedIndices > &indices)
 Performs static condensation on a square matrix. More...
 
template<typename Derived , size_t sizeOfRemovedCols>
auto Ikarus::removeCol (const Eigen::MatrixBase< Derived > &E, const std::array< size_t, sizeOfRemovedCols > &indices)
 Removes specified columns from a matrix. More...
 
template<typename ST , typename MaterialImpl >
auto Ikarus::toVoigtAndMaybeReduce (const Eigen::Matrix< ST, 3, 3 > &E, const MaterialImpl &material, bool isStrain=true)
 Converts a 3x3 matrix to Voigt notation, possibly reducing it based on material properties. More...
 
template<typename M , typename Derived >
decltype(auto) Ikarus::enlargeIfReduced (const Eigen::MatrixBase< Derived > &E)
 Enlarges a matrix if it reduced in the context of material laws, i.e., VanishingStress If the material is not reduced the untouched matrix is returned and rendering the function as a NoOp. More...
 
template<typename MessageType >
MessageType & Ikarus::increment (MessageType &e)
 Increments the given enum value. More...
 
std::tuple< Dune::Functions::Polynomial< double >, double > Ikarus::utils::polyfit (const Eigen::Ref< const Eigen::VectorXd > &x, const Eigen::Ref< const Eigen::VectorXd > &y, int deg)
 Fits a polynomial of a given degree to the given data points. More...
 

Variables

template<int dim>
constexpr auto Ikarus::voigtNotationContainer = std::get<dim - 1>(Impl::voigtIndices)
 Container for Voigt notation indices based on dimension.1D: 0,0 2D: 0,0; 1,1; 0,1 3D: 0,0; 1,1; 2,2; 1,2; 0,2; 0,1. More...
 

Detailed Description

Macro Definition Documentation

◆ MAKE_ENUM

#define MAKE_ENUM (   type,
  ... 
)
Value:
enum class type \
{ \
BEGIN, \
__VA_ARGS__, \
END \
}; \
constexpr std::string toString(type _e) { \
using enum type; \
switch (_e) { \
ENUM_CASE(BEGIN) \
FOR_EACH(ENUM_CASE, __VA_ARGS__) \
ENUM_CASE(END) \
} \
__builtin_unreachable(); \
}
#define ENUM_CASE(name)
Definition: makeenum.hh:27
constexpr std::string toString(ScalarAffordances _e)
Definition: ferequirements.hh:37
Parameters
typeThe name of the enum class.
...Enumerators to be included between BEGIN and END.

Function Documentation

◆ addInEmbedding()

template<typename ManifoldPoint , typename Derived >
requires (Ikarus::Concepts::AddAssignAble<ManifoldPoint, decltype(b.template segment<ManifoldPoint::valueSize>(0))> and requires() { ManifoldPoint::valueSize; })
Dune::BlockVector< ManifoldPoint > & Ikarus::addInEmbedding ( Dune::BlockVector< ManifoldPoint > &  a,
const Eigen::MatrixBase< Derived > &  b 
)
Template Parameters
ManifoldPointManifold type.
DerivedManifoldPoint of the input Eigen matrix.
Parameters
aInput Dune::BlockVector.
bInput Eigen::Matrix.
Returns
Reference to the modified Dune::BlockVector.

◆ checkGradient()

template<typename NonlinearOperator , typename UpdateType = typename NonlinearOperator::template ParameterValue<0>>
bool Ikarus::utils::checkGradient ( NonlinearOperator &  nonLinOp,
CheckFlags  checkFlags = CheckFlags(),
std::function< void(typename NonlinearOperator::template ParameterValue< 0 > &, const UpdateType &)>  p_updateFunction = [](typename NonlinearOperator::template ParameterValue<0>& a, const UpdateType& b) { a += b; } 
)

The checkgradient function is inspired by http://sma.epfl.ch/~nboumal/book/ Chapter 4.8 and https://github.com/NicolasBoumal/manopt/blob/master/manopt/tools/checkdiff.m

Template Parameters
NonlinearOperatorType of the nonlinear operator.
UpdateTypeType of the update.
Parameters
nonLinOpThe nonlinear operator.
checkFlagsFlags for the check.
p_updateFunctionUpdate function.
Returns
True if the check passed, false otherwise.

◆ checkHessian()

template<typename NonlinearOperator , typename UpdateType = typename NonlinearOperator::template ParameterValue<0>>
bool Ikarus::utils::checkHessian ( NonlinearOperator &  nonLinOp,
CheckFlags  checkFlags = CheckFlags(),
std::function< void(typename NonlinearOperator::template ParameterValue< 0 > &, const UpdateType &)>  p_updateFunction = [](typename NonlinearOperator::template ParameterValue<0>& a, const UpdateType& b) { a += b; } 
)

The checkHessian function is inspired by http://sma.epfl.ch/~nboumal/book/ Chapter 6.8 and https://github.com/NicolasBoumal/manopt/blob/master/manopt/tools/checkhessian.m

Template Parameters
NonlinearOperatorType of the nonlinear operator.
UpdateTypeType of the update.
Parameters
nonLinOpThe nonlinear operator.
checkFlagsFlags for the check.
p_updateFunctionUpdate function.
Returns
True if the check passed, false otherwise.

◆ checkJacobian()

template<typename NonlinearOperator , typename UpdateType = typename NonlinearOperator::template ParameterValue<0>>
bool Ikarus::utils::checkJacobian ( NonlinearOperator &  nonLinOp,
CheckFlags  checkFlags = CheckFlags(),
std::function< void(typename NonlinearOperator::template ParameterValue< 0 > &, const UpdateType &)>  p_updateFunction = [](typename NonlinearOperator::template ParameterValue<0>& a, const UpdateType& b) { a += b; } 
)

The checkjacobian function is inspired by http://sma.epfl.ch/~nboumal/book/ Chapter 4.8 and https://github.com/NicolasBoumal/manopt/blob/master/manopt/tools/checkdiff.m

Template Parameters
NonlinearOperatorType of the nonlinear operator.
UpdateTypeType of the update.
Parameters
nonLinOpThe nonlinear operator.
checkFlagsFlags for the check.
p_updateFunctionUpdate function.
Returns
True if the check passed, false otherwise.

◆ correctionSize()

template<typename Type >
requires requires { Type::correctionSize; }
size_t Ikarus::correctionSize ( const Dune::BlockVector< Type > &  a)
Template Parameters
TypeManifold type.
Parameters
aInput Dune::BlockVector.
Returns
Total correction size.
Here is the caller graph for this function:

◆ createRandomVector()

template<typename FieldVectorT >
auto Ikarus::createRandomVector ( typename FieldVectorT::value_type  lower = -1,
typename FieldVectorT::value_type  upper = 1 
)
Template Parameters
FieldVectorTThe type of the vector.
Parameters
lowerThe lower bound of the random values (default is -1).
upperThe upper bound of the random values (default is 1).
Returns
A random vector within the specified range.

◆ enlargeIfReduced()

template<typename M , typename Derived >
decltype(auto) Ikarus::enlargeIfReduced ( const Eigen::MatrixBase< Derived > &  E)
Template Parameters
MType of the material.
DerivedType of the input matrix.
Parameters
EInput matrix.
Returns
auto Resulting matrix based on material properties.

This function takes an input matrix and, based on the material properties, either returns the original matrix (if it is not reduced) or enlarges the matrix by filling in the specified columns with zeros (if it is reduced).

◆ findLineSegment()

std::tuple< Dune::Functions::Polynomial< double >, decltype(Eigen::seq(0, 0))> Ikarus::utils::findLineSegment ( const Eigen::VectorXd &  x,
const Eigen::VectorXd &  y,
int  segmentSize 
)

his function is inspired by the MATLAB code at: https://github.com/NicolasBoumal/manopt/blob/master/manopt/tools/identify_linear_piece.m It is designed to find the most linear segment in a set of data points

Parameters
xThe x-coordinates of the data points.
yThe y-coordinates of the data points.
segmentSizeThe size of the line segment to be identified.
Returns
A tuple containing the polynomial representing the identified line segment and the indices of the data points in the segment.

◆ flatPreBasis()

template<class PreBasis >
decltype(auto) Ikarus::flatPreBasis ( const PreBasis &  preBasis)

◆ hessianN()

template<typename Fun , typename... Vars, typename... Args, typename U , typename G , typename H >
void Ikarus::utils::hessianN ( const Fun &  f,
const autodiff::Wrt< Vars... > &  wrt,
const autodiff::At< Args... > &  at,
U &  u,
std::array< G, U::RowsAtCompileTime > &  g,
std::array< H, U::RowsAtCompileTime > &  h 
)
Template Parameters
FunThe type of the function to be differentiated.
VarsThe types of the variables with respect to which the Hessian is computed.
ArgsThe types of the arguments passed to the function.
UThe type representing the result of the function evaluation.
GThe type representing the gradient of the function.
HThe type representing the Hessian matrix.
Parameters
fThe function to be differentiated.
wrtThe variables with respect to which the Hessian is computed.
atThe values at which the Hessian is evaluated.
uThe result of the function evaluation.
gThe gradient of the function.
hThe Hessian matrix (output).

◆ increment()

template<typename MessageType >
MessageType & Ikarus::increment ( MessageType &  e)
Template Parameters
MessageTypeThe enum type.
Parameters
eThe enum value to increment.
Returns
The incremented enum value.
Exceptions
Dune::RangeErrorif trying to increment MessageType::END.
Here is the caller graph for this function:

◆ norm() [1/2]

template<typename Derived >
requires (!std::floating_point<Derived>)
auto Ikarus::norm ( const Eigen::MatrixBase< Derived > &  v)
Template Parameters
DerivedType of the input Eigen matrix.
Parameters
vInput Eigen matrix.
Returns
Norm of the matrix.

◆ norm() [2/2]

auto Ikarus::norm ( const std::floating_point auto &  v)
Parameters
vInput scalar.
Returns
Absolute value of the scalar.

◆ obtainLagrangeNodePositions()

template<int size, typename LV >
void Ikarus::utils::obtainLagrangeNodePositions ( const LV &  localView,
std::vector< Dune::FieldVector< double, size > > &  lagrangeNodeCoords 
)
Template Parameters
sizeSize of the nodal coordinate vector
LVType of the local view
Parameters
localViewLocal view bounded to an element
lagrangeNodeCoordsA vector of nodal coordinates to be updated

◆ operator*()

template<typename Scalar , int size>
auto Ikarus::operator* ( const Eigen::DiagonalMatrix< Scalar, size > &  a,
const Eigen::DiagonalMatrix< Scalar, size > &  b 
)
Template Parameters
ScalarScalar type.
sizeSize of the diagonal matrix.
Parameters
aInput DiagonalMatrix.
bInput DiagonalMatrix.
Returns
Product of the two DiagonalMatrices.

◆ operator+() [1/4]

template<typename Derived , typename Scalar , int size>
auto Ikarus::operator+ ( const Eigen::DiagonalMatrix< Scalar, size > &  a,
const Eigen::MatrixBase< Derived > &  b 
)
Template Parameters
DerivedType of the input Eigen matrix.
ScalarScalar type.
sizeSize of the diagonal matrix.
Parameters
aInput DiagonalMatrix.
bInput Eigen matrix.
Returns
Sum of the DiagonalMatrix and Eigen matrix.

◆ operator+() [2/4]

template<typename Derived , typename Derived2 >
auto Ikarus::operator+ ( const Eigen::DiagonalWrapper< Derived > &  a,
const Eigen::MatrixBase< Derived2 > &  b 
)
Template Parameters
DerivedType of the input Eigen DiagonalWrapper.
Derived2Type of the input Eigen matrix.
Parameters
aInput Eigen DiagonalWrapper.
bInput Eigen matrix.
Returns
Sum of DiagonalWrapper and Eigen matrix.

◆ operator+() [3/4]

template<typename Derived , typename Scalar , int size>
auto Ikarus::operator+ ( const Eigen::MatrixBase< Derived > &  a,
const Eigen::DiagonalMatrix< Scalar, size > &  b 
)
Template Parameters
DerivedType of the input Eigen matrix.
ScalarScalar type.
sizeSize of the diagonal matrix.
Parameters
aInput Eigen matrix.
bInput DiagonalMatrix.
Returns
Sum of the Eigen matrix and DiagonalMatrix.

◆ operator+() [4/4]

template<typename Derived , typename Derived2 >
auto Ikarus::operator+ ( const Eigen::MatrixBase< Derived > &  a,
const Eigen::DiagonalWrapper< Derived2 > &  b 
)
Template Parameters
DerivedType of the input Eigen matrix.
Derived2Type of the input Eigen DiagonalWrapper.
Parameters
aInput Eigen matrix.
bInput Eigen DiagonalWrapper.
Returns
Sum of Eigen matrix and DiagonalWrapper.

◆ operator+=() [1/3]

template<typename Type , typename Derived >
requires (Ikarus::Concepts::AddAssignAble<Type, decltype(b.template segment<Type::correctionSize>(0))> and requires() { Type::correctionSize; })
Dune::BlockVector< Type > & Ikarus::operator+= ( Dune::BlockVector< Type > &  a,
const Eigen::MatrixBase< Derived > &  b 
)
Template Parameters
TypeManifold type.
DerivedType of the input Eigen matrix.
Parameters
aInput Dune::BlockVector.
bInput Eigen::Matrix.
Returns
Reference to the modified Dune::BlockVector.

◆ operator+=() [2/3]

template<typename... Types, typename Derived >
Dune::TupleVector< Types... > & Ikarus::operator+= ( Dune::TupleVector< Types... > &  a,
const Eigen::MatrixBase< Derived > &  b 
)
Template Parameters
TypesTypes of the elements in the TupleVector.
DerivedType of the input Eigen matrix.
Parameters
aInput Dune::TupleVector.
bInput Eigen::Matrix.
Returns
Reference to the modified Dune::TupleVector.

◆ operator+=() [3/3]

template<typename Scalar , int size>
auto Ikarus::operator+= ( Eigen::DiagonalMatrix< Scalar, size > &  a,
const Eigen::DiagonalMatrix< Scalar, size > &  b 
)
Template Parameters
ScalarScalar type.
sizeSize of the diagonal matrix.
Parameters
aInput DiagonalMatrix.
bInput DiagonalMatrix.
Returns
Reference to the modified DiagonalMatrix.

◆ operator-()

template<typename Scalar , int size>
auto Ikarus::operator- ( const Eigen::DiagonalMatrix< Scalar, size > &  a)
Template Parameters
ScalarScalar type.
sizeSize of the diagonal matrix.
Parameters
aInput DiagonalMatrix.
Returns
Negation of the DiagonalMatrix.

◆ operator-=()

template<typename Type , typename Derived >
requires (Ikarus::Concepts::AddAssignAble<Type, decltype(b.template segment<Type::correctionSize>(0))> and requires() { Type::correctionSize; })
Dune::BlockVector< Type > & Ikarus::operator-= ( Dune::BlockVector< Type > &  a,
const Eigen::MatrixBase< Derived > &  b 
)
Template Parameters
TypeManifold type.
DerivedType of the input Eigen matrix.
Parameters
aInput Dune::BlockVector.
bInput Eigen::Matrix.
Returns
Reference to the modified Dune::BlockVector.

◆ operator<<()

template<typename Scalar , int size>
std::ostream & Ikarus::operator<< ( std::ostream &  os,
const Eigen::DiagonalMatrix< Scalar, size > &  a 
)
Template Parameters
ScalarScalar type.
sizeSize of the diagonal matrix.
Parameters
osOutput stream.
aInput DiagonalMatrix.
Returns
Reference to the output stream.

◆ orthonormalizeMatrixColumns()

template<typename Derived >
auto Ikarus::orthonormalizeMatrixColumns ( const Eigen::MatrixBase< Derived > &  A)
Template Parameters
DerivedType of the input matrix.
Parameters
AThe input matrix.
Returns
Eigen Matrix with orthonormalized columns.

◆ polyfit()

std::tuple< Dune::Functions::Polynomial< double >, double > Ikarus::utils::polyfit ( const Eigen::Ref< const Eigen::VectorXd > &  x,
const Eigen::Ref< const Eigen::VectorXd > &  y,
int  deg 
)
Parameters
xThe input vector of x-coordinates.
yThe input vector of y-coordinates.
degThe degree of the polynomial to fit.
Returns
std::tuple<Dune::Functions::Polynomial<double>, double> A tuple containing the fitted polynomial and the least square error.

◆ printForMaple()

template<typename Derived >
void Ikarus::printForMaple ( const Eigen::EigenBase< Derived > &  A)
Template Parameters
DerivedThe derived type of the matrix.
Parameters
AThe input matrix.

◆ referenceElementSubEntityPositions()

template<typename FE >
auto Ikarus::utils::referenceElementSubEntityPositions ( FE fe,
int  codim 
)
Template Parameters
FEType of the finite element
Parameters
fefinite element
codimcodim of requested subentity
Returns
view over the position of the subenties
Here is the caller graph for this function:

◆ referenceElementVertexPositions()

template<typename FE >
auto Ikarus::utils::referenceElementVertexPositions ( FE fe)
Template Parameters
FEType of the finite element
Parameters
fefinite element
Returns

◆ removeCol()

template<typename Derived , size_t sizeOfRemovedCols>
auto Ikarus::removeCol ( const Eigen::MatrixBase< Derived > &  E,
const std::array< size_t, sizeOfRemovedCols > &  indices 
)
Template Parameters
DerivedType of the matrix.
sizeOfRemovedColsSize of the columns to be removed.
Parameters
EInput matrix.
indicesArray of column indices to be removed.
Returns
Resulting matrix after removing specified columns.

This function removes specified columns from a matrix. It computes the remaining columns after removing the specified indices and returns the resulting matrix.

Here is the caller graph for this function:

◆ skew() [1/2]

template<typename Derived >
Derived Ikarus::skew ( const Eigen::MatrixBase< Derived > &  A)
Template Parameters
DerivedType of the input Eigen matrix.
Parameters
AInput Eigen matrix.
Returns
Skew part of the matrix.
Here is the caller graph for this function:

◆ skew() [2/2]

template<typename ScalarType >
Eigen::Matrix< ScalarType, 3, 3 > Ikarus::skew ( const Eigen::Vector< ScalarType, 3 > &  a)
Template Parameters
ScalarTypeThe type of the coordinates in the vector.
Parameters
aThe vector.
Returns
The skew matrix.

◆ staticCondensation()

template<typename Derived , size_t sizeOfCondensedIndices>
auto Ikarus::staticCondensation ( const Eigen::MatrixBase< Derived > &  E,
const std::array< size_t, sizeOfCondensedIndices > &  indices 
)
Template Parameters
DerivedType of the matrix.
sizeOfCondensedIndicesSize of the condensed indices.
Parameters
EInput matrix.
indicesArray of indices to be condensed.
Returns
Resulting matrix after static condensation.

This function performs static condensation on a square matrix. It removes the specified indices from the matrix, computes the remaining submatrices (K11, K12, K22), and returns the result of the static condensation.

Here is the caller graph for this function:

◆ sym()

template<typename Derived >
Derived Ikarus::sym ( const Eigen::MatrixBase< Derived > &  A)
Template Parameters
DerivedType of the input Eigen matrix.
Parameters
AInput Eigen matrix.
Returns
Symmetric part of the matrix.

◆ toVoigtAndMaybeReduce()

template<typename ST , typename MaterialImpl >
auto Ikarus::toVoigtAndMaybeReduce ( const Eigen::Matrix< ST, 3, 3 > &  E,
const MaterialImpl &  material,
bool  isStrain = true 
)
Template Parameters
STScalar type of the matrix.
MaterialImplType of the material implementation.
Parameters
EInput 3x3 matrix.
materialReference to the material implementation.
isStrainFlag indicating if the matrix represents strain (default is true).
Returns
Resulting matrix in Voigt notation.

This function converts a 3x3 matrix to its Voigt notation. If the material is not reduced, the full Voigt notation is returned. Otherwise, the specified columns (based on material properties, such as VanishingStress) are removed, and the reduced Voigt notation is returned.

◆ valueSize()

template<typename Type >
requires requires { Type::valueSize; }
size_t Ikarus::valueSize ( const Dune::BlockVector< Type > &  a)
Template Parameters
TypeManifold type.
Parameters
aInput Dune::BlockVector.
Returns
Total value size.
Here is the caller graph for this function:

◆ viewAsEigenMatrixAsDynFixed() [1/2]

template<typename ValueType >
auto Ikarus::viewAsEigenMatrixAsDynFixed ( const Dune::BlockVector< ValueType > &  blockedVector)
Template Parameters
ValueTypeType of the elements in the BlockVector.
Parameters
blockedVectorInput Dune::BlockVector.
Returns
Eigen::Map of the BlockVector as a dynamic Eigen::Matrix (const version).

◆ viewAsEigenMatrixAsDynFixed() [2/2]

template<typename ValueType >
auto Ikarus::viewAsEigenMatrixAsDynFixed ( Dune::BlockVector< ValueType > &  blockedVector)
Template Parameters
ValueTypeType of the elements in the BlockVector.
Parameters
blockedVectorInput Dune::BlockVector.
Returns
Eigen::Map of the BlockVector as a dynamic Eigen::Matrix.
Here is the caller graph for this function:

◆ viewAsEigenMatrixFixedDyn() [1/2]

template<typename ValueType >
auto Ikarus::viewAsEigenMatrixFixedDyn ( const Dune::BlockVector< ValueType > &  blockedVector)
Template Parameters
ValueTypeType of the elements in the BlockVector.
Parameters
blockedVectorInput Dune::BlockVector.
Returns
Eigen::Map of the BlockVector as a fixed-size Eigen::Matrix with dynamic columns (const version).

◆ viewAsEigenMatrixFixedDyn() [2/2]

template<typename ValueType >
auto Ikarus::viewAsEigenMatrixFixedDyn ( Dune::BlockVector< ValueType > &  blockedVector)
Template Parameters
ValueTypeType of the elements in the BlockVector.
Parameters
blockedVectorInput Dune::BlockVector.
Returns
Eigen::Map of the BlockVector as a fixed-size Eigen::Matrix with dynamic columns.

◆ viewAsFlatEigenVector() [1/2]

template<typename ValueType >
auto Ikarus::viewAsFlatEigenVector ( const Dune::BlockVector< ValueType > &  blockedVector)
Template Parameters
ValueTypeType of the elements in the BlockVector.
Parameters
blockedVectorInput Dune::BlockVector.
Returns
Eigen::Map of the BlockVector as a flat Eigen::Vector (const version).

◆ viewAsFlatEigenVector() [2/2]

template<typename ValueType >
auto Ikarus::viewAsFlatEigenVector ( Dune::BlockVector< ValueType > &  blockedVector)
Template Parameters
ValueTypeType of the elements in the BlockVector.
Parameters
blockedVectorInput Dune::BlockVector.
Returns
Eigen::Map of the BlockVector as a flat Eigen::Vector.
Here is the caller graph for this function:

Variable Documentation

◆ voigtNotationContainer

template<int dim>
constexpr auto Ikarus::voigtNotationContainer = std::get<dim - 1>(Impl::voigtIndices)
constexpr
Template Parameters
dimThe dimension for which Voigt indices are needed.