Collection of several utilities. More...
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::obtainLagrangeGlobalNodePositions (const LV &localView, std::vector< Dune::FieldVector< double, size > > &lagrangeNodeGlobalCoords) |
A function to obtain the global positions of the nodes of an element with Lagrangian basis. More... | |
template<int size, typename Basis > | |
auto | Ikarus::utils::globalIndexFromGlobalPosition (const Basis &basis, const Dune::FieldVector< double, size > &pos) |
A helper function to obtain the global index from the global positions for a Lagrange node. 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... | |
void | addBindingsToUtils () |
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... | |
#define MAKE_ENUM | ( | type, | |
... | |||
) |
type | The name of the enum class. |
... | Enumerators to be included between BEGIN and END. |
void addBindingsToUtils | ( | ) |
Converts a square 1x1, 2x2 or 3x3 matrix to a Voigt notation vector.
E | Input matrix of size (size x size). |
isStrain | Flag indicating whether the conversion is for strain (true) or not (false) (default is true).. |
This function converts a square matrix to a Voigt notation vector, which contains the unique components of the input matrix.
The optional isStrain parameter allows the user to specify whether the conversion is intended for strain calculations. If isStrain is true, the off-diagonal components are multiplied by 2, providing the correct Voigt notation for symmetric strain tensors.
Converts a vector given in Voigt notation to a matrix.
EVoigt | Voigt notation vector. |
isStrain | Flag indicating whether the vector represents a strain (default is true). |
This function converts a vector given in Voigt notation to the corresponding matrix. The conversion depends on the size The parameter isStrain
is used to determine the conversion factor for off-diagonal components, which need to be divided by 2 in the matrix representation if the quantity is a strain tensor.
The function requires that the size of the Voigt notation vector is valid (1, 3, or 6).
Dune::BlockVector< ManifoldPoint > & Ikarus::addInEmbedding | ( | Dune::BlockVector< ManifoldPoint > & | a, |
const Eigen::MatrixBase< Derived > & | b | ||
) |
ManifoldPoint | Manifold type. |
Derived | ManifoldPoint of the input Eigen matrix. |
a | Input Dune::BlockVector. |
b | Input Eigen::Matrix. |
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
NonlinearOperator | Type of the nonlinear operator. |
UpdateType | Type of the update. |
nonLinOp | The nonlinear operator. |
checkFlags | Flags for the check. |
p_updateFunction | Update function. |
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
NonlinearOperator | Type of the nonlinear operator. |
UpdateType | Type of the update. |
nonLinOp | The nonlinear operator. |
checkFlags | Flags for the check. |
p_updateFunction | Update function. |
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
NonlinearOperator | Type of the nonlinear operator. |
UpdateType | Type of the update. |
nonLinOp | The nonlinear operator. |
checkFlags | Flags for the check. |
p_updateFunction | Update function. |
size_t Ikarus::correctionSize | ( | const Dune::BlockVector< Type > & | a | ) |
Type | Manifold type. |
a | Input Dune::BlockVector. |
auto Ikarus::createRandomVector | ( | typename FieldVectorT::value_type | lower = -1 , |
typename FieldVectorT::value_type | upper = 1 |
||
) |
FieldVectorT | The type of the vector. |
lower | The lower bound of the random values (default is -1). |
upper | The upper bound of the random values (default is 1). |
decltype(auto) Ikarus::enlargeIfReduced | ( | const Eigen::MatrixBase< Derived > & | E | ) |
M | Type of the material. |
Derived | Type of the input matrix. |
E | Input matrix. |
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).
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
x | The x-coordinates of the data points. |
y | The y-coordinates of the data points. |
segmentSize | The size of the line segment to be identified. |
decltype(auto) Ikarus::flatPreBasis | ( | const PreBasis & | preBasis | ) |
auto Ikarus::utils::globalIndexFromGlobalPosition | ( | const Basis & | basis, |
const Dune::FieldVector< double, size > & | pos | ||
) |
size | Size of the global nodal coordinate vector |
Basis | Type of the basis. |
basis | The grid basis. |
pos | Global position |
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 | ||
) |
Fun | The type of the function to be differentiated. |
Vars | The types of the variables with respect to which the Hessian is computed. |
Args | The types of the arguments passed to the function. |
U | The type representing the result of the function evaluation. |
G | The type representing the gradient of the function. |
H | The type representing the Hessian matrix. |
f | The function to be differentiated. |
wrt | The variables with respect to which the Hessian is computed. |
at | The values at which the Hessian is evaluated. |
u | The result of the function evaluation. |
g | The gradient of the function. |
h | The Hessian matrix (output). |
MessageType & Ikarus::increment | ( | MessageType & | e | ) |
MessageType | The enum type. |
e | The enum value to increment. |
Dune::RangeError | if trying to increment MessageType::END. |
auto Ikarus::norm | ( | const Eigen::MatrixBase< Derived > & | v | ) |
auto Ikarus::norm | ( | const std::floating_point auto & | v | ) |
v | Input scalar. |
void Ikarus::utils::obtainLagrangeGlobalNodePositions | ( | const LV & | localView, |
std::vector< Dune::FieldVector< double, size > > & | lagrangeNodeGlobalCoords | ||
) |
size | Size of the global nodal coordinate vector |
LV | Type of the local view |
localView | Local view bounded to an element |
lagrangeNodeGlobalCoords | A vector of global nodal coordinates to be updated |
auto Ikarus::operator* | ( | const Eigen::DiagonalMatrix< Scalar, size > & | a, |
const Eigen::DiagonalMatrix< Scalar, size > & | b | ||
) |
Scalar | Scalar type. |
size | Size of the diagonal matrix. |
a | Input DiagonalMatrix. |
b | Input DiagonalMatrix. |
auto Ikarus::operator+ | ( | const Eigen::DiagonalMatrix< Scalar, size > & | a, |
const Eigen::MatrixBase< Derived > & | b | ||
) |
auto Ikarus::operator+ | ( | const Eigen::DiagonalWrapper< Derived > & | a, |
const Eigen::MatrixBase< Derived2 > & | b | ||
) |
auto Ikarus::operator+ | ( | const Eigen::MatrixBase< Derived > & | a, |
const Eigen::DiagonalMatrix< Scalar, size > & | b | ||
) |
auto Ikarus::operator+ | ( | const Eigen::MatrixBase< Derived > & | a, |
const Eigen::DiagonalWrapper< Derived2 > & | b | ||
) |
Dune::BlockVector< Type > & Ikarus::operator+= | ( | Dune::BlockVector< Type > & | a, |
const Eigen::MatrixBase< Derived > & | b | ||
) |
Type | Manifold type. |
Derived | Type of the input Eigen matrix. |
a | Input Dune::BlockVector. |
b | Input Eigen::Matrix. |
Dune::TupleVector< Types... > & Ikarus::operator+= | ( | Dune::TupleVector< Types... > & | a, |
const Eigen::MatrixBase< Derived > & | b | ||
) |
Types | Types of the elements in the TupleVector. |
Derived | Type of the input Eigen matrix. |
a | Input Dune::TupleVector. |
b | Input Eigen::Matrix. |
auto Ikarus::operator+= | ( | Eigen::DiagonalMatrix< Scalar, size > & | a, |
const Eigen::DiagonalMatrix< Scalar, size > & | b | ||
) |
Scalar | Scalar type. |
size | Size of the diagonal matrix. |
a | Input DiagonalMatrix. |
b | Input DiagonalMatrix. |
auto Ikarus::operator- | ( | const Eigen::DiagonalMatrix< Scalar, size > & | a | ) |
Scalar | Scalar type. |
size | Size of the diagonal matrix. |
a | Input DiagonalMatrix. |
Dune::BlockVector< Type > & Ikarus::operator-= | ( | Dune::BlockVector< Type > & | a, |
const Eigen::MatrixBase< Derived > & | b | ||
) |
Type | Manifold type. |
Derived | Type of the input Eigen matrix. |
a | Input Dune::BlockVector. |
b | Input Eigen::Matrix. |
std::ostream & Ikarus::operator<< | ( | std::ostream & | os, |
const Eigen::DiagonalMatrix< Scalar, size > & | a | ||
) |
Scalar | Scalar type. |
size | Size of the diagonal matrix. |
os | Output stream. |
a | Input DiagonalMatrix. |
auto Ikarus::orthonormalizeMatrixColumns | ( | const Eigen::MatrixBase< Derived > & | A | ) |
Derived | Type of the input matrix. |
A | The input matrix. |
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 | ||
) |
x | The input vector of x-coordinates. |
y | The input vector of y-coordinates. |
deg | The degree of the polynomial to fit. |
void Ikarus::printForMaple | ( | const Eigen::EigenBase< Derived > & | A | ) |
Derived | The derived type of the matrix. |
A | The input matrix. |
auto Ikarus::utils::referenceElementSubEntityPositions | ( | FE & | fe, |
int | codim | ||
) |
FE | Type of the finite element |
fe | finite element |
codim | codim of requested subentity |
auto Ikarus::removeCol | ( | const Eigen::MatrixBase< Derived > & | E, |
const std::array< size_t, sizeOfRemovedCols > & | indices | ||
) |
Derived | Type of the matrix. |
sizeOfRemovedCols | Size of the columns to be removed. |
E | Input matrix. |
indices | Array of column indices to be removed. |
This function removes specified columns from a matrix. It computes the remaining columns after removing the specified indices and returns the resulting matrix.
Derived Ikarus::skew | ( | const Eigen::MatrixBase< Derived > & | A | ) |
Eigen::Matrix< ScalarType, 3, 3 > Ikarus::skew | ( | const Eigen::Vector< ScalarType, 3 > & | a | ) |
ScalarType | The type of the coordinates in the vector. |
a | The vector. |
auto Ikarus::staticCondensation | ( | const Eigen::MatrixBase< Derived > & | E, |
const std::array< size_t, sizeOfCondensedIndices > & | indices | ||
) |
Derived | Type of the matrix. |
sizeOfCondensedIndices | Size of the condensed indices. |
E | Input matrix. |
indices | Array of indices to be condensed. |
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.
Derived Ikarus::sym | ( | const Eigen::MatrixBase< Derived > & | A | ) |
auto Ikarus::toVoigtAndMaybeReduce | ( | const Eigen::Matrix< ST, 3, 3 > & | E, |
const MaterialImpl & | material, | ||
bool | isStrain = true |
||
) |
ST | Scalar type of the matrix. |
MaterialImpl | Type of the material implementation. |
E | Input 3x3 matrix. |
material | Reference to the material implementation. |
isStrain | Flag indicating if the matrix represents strain (default is true). |
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.
size_t Ikarus::valueSize | ( | const Dune::BlockVector< Type > & | a | ) |
Type | Manifold type. |
a | Input Dune::BlockVector. |
auto Ikarus::viewAsEigenMatrixAsDynFixed | ( | const Dune::BlockVector< ValueType > & | blockedVector | ) |
ValueType | Type of the elements in the BlockVector. |
blockedVector | Input Dune::BlockVector. |
auto Ikarus::viewAsEigenMatrixAsDynFixed | ( | Dune::BlockVector< ValueType > & | blockedVector | ) |
ValueType | Type of the elements in the BlockVector. |
blockedVector | Input Dune::BlockVector. |
auto Ikarus::viewAsEigenMatrixFixedDyn | ( | const Dune::BlockVector< ValueType > & | blockedVector | ) |
ValueType | Type of the elements in the BlockVector. |
blockedVector | Input Dune::BlockVector. |
auto Ikarus::viewAsEigenMatrixFixedDyn | ( | Dune::BlockVector< ValueType > & | blockedVector | ) |
ValueType | Type of the elements in the BlockVector. |
blockedVector | Input Dune::BlockVector. |
auto Ikarus::viewAsFlatEigenVector | ( | const Dune::BlockVector< ValueType > & | blockedVector | ) |
ValueType | Type of the elements in the BlockVector. |
blockedVector | Input Dune::BlockVector. |
auto Ikarus::viewAsFlatEigenVector | ( | Dune::BlockVector< ValueType > & | blockedVector | ) |
ValueType | Type of the elements in the BlockVector. |
blockedVector | Input Dune::BlockVector. |
|
constexpr |
dim | The dimension for which Voigt indices are needed. |