KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
List of all members
Kratos::MathUtils< TDataType > Class Template Reference

Various mathematical utilitiy functions. More...

#include <math_utils.h>

Collaboration diagram for Kratos::MathUtils< TDataType >:

Type Definitions

using MatrixType = Matrix
 The matrix type. More...
 
using VectorType = Vector
 The vector type. More...
 
using SizeType = std::size_t
 The size type. More...
 
using IndexType = std::size_t
 The index type. More...
 
using IndirectArrayType = boost::numeric::ublas::indirect_array< DenseVector< std::size_t > >
 The indirect array type. More...
 
static constexpr double ZeroTolerance = std::numeric_limits<double>::epsilon()
 The machine precision. More...
 

Operations

 TDim
 
TDimrA
 
TDim BoundedMatrix< double, TDim, TDim > & rEigenVectorsMatrix
 
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > & rEigenValuesMatrix
 
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > const double Tolerance = 1.0e-18
 
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > const double const SizeType MaxIterations
 
static double GetZeroTolerance ()
 This function returns the machine precision. More...
 
template<bool TCheck>
static double Heron (double a, double b, double c)
 In geometry, Heron's formula (sometimes called Hero's formula), named after Hero of Alexandria, gives the area of a triangle by requiring no arbitrary choice of side as base or vertex as origin, contrary to other formulas for the area of a triangle, such as half the base times the height or half the norm of a cross product of two sides. More...
 
template<class TMatrixType >
static double Cofactor (const TMatrixType &rMat, IndexType i, IndexType j)
 Calculates the cofactor. More...
 
template<class TMatrixType >
static MatrixType CofactorMatrix (const TMatrixType &rMat)
 Calculates the cofactor matrix. More...
 
template<SizeType TDim>
static BoundedMatrix< double, TDim, TDimInvertMatrix (const BoundedMatrix< double, TDim, TDim > &rInputMatrix, double &rInputMatrixDet, const double Tolerance=ZeroTolerance)
 Calculates the inverse of a 2x2, 3x3 and 4x4 matrices (using bounded matrix for performance) More...
 
template<class TMatrix1 , class TMatrix2 >
static bool CheckConditionNumber (const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, const double Tolerance=std::numeric_limits< double >::epsilon(), const bool ThrowError=true)
 This method checks the condition number of amtrix. More...
 
template<class TMatrix1 , class TMatrix2 >
static void GeneralizedInvertMatrix (const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, double &rInputMatrixDet, const double Tolerance=ZeroTolerance)
 It inverts non square matrices (https://en.wikipedia.org/wiki/Inverse_element#Matrices) More...
 
static void Solve (MatrixType A, VectorType &rX, const VectorType &rB)
 This function is designed to be called when a dense linear system is needed to be solved. More...
 
template<class TMatrix1 , class TMatrix2 >
static void InvertMatrix (const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, double &rInputMatrixDet, const double Tolerance=ZeroTolerance)
 It inverts matrices of order 2, 3 and 4. More...
 
template<class TMatrix1 , class TMatrix2 >
static void InvertMatrix2 (const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, double &rInputMatrixDet)
 It inverts matrices of order 2. More...
 
template<class TMatrix1 , class TMatrix2 >
static void InvertMatrix3 (const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, double &rInputMatrixDet)
 It inverts matrices of order 3. More...
 
template<class TMatrix1 , class TMatrix2 >
static void InvertMatrix4 (const TMatrix1 &rInputMatrix, TMatrix2 &rInvertedMatrix, double &rInputMatrixDet)
 It inverts matrices of order 4. More...
 
template<class TMatrixType >
static double Det2 (const TMatrixType &rA)
 Calculates the determinant of a matrix of dimension 2x2 (no check performed on matrix size) More...
 
template<class TMatrixType >
static double Det3 (const TMatrixType &rA)
 Calculates the determinant of a matrix of dimension 3*3 (no check performed on matrix size) More...
 
template<class TMatrixType >
static double Det4 (const TMatrixType &rA)
 Calculates the determinant of a matrix of dimension 4*4 (no check performed on matrix size) More...
 
template<class TMatrixType >
static double Det (const TMatrixType &rA)
 Calculates the determinant of a matrix of a square matrix of any size (no check performed on release mode) More...
 
template<class TMatrixType >
static double GeneralizedDet (const TMatrixType &rA)
 Calculates the determinant of any matrix (no check performed on matrix size) More...
 
static double Dot3 (const Vector &a, const Vector &b)
 Performs the dot product of two vectors of dimension 3. More...
 
static double Dot (const Vector &rFirstVector, const Vector &rSecondVector)
 Performs the dot product of two vectors of arbitrary size. More...
 
template<class TVectorType >
static double Norm3 (const TVectorType &a)
 Calculates the norm of vector "a" which is assumed to be of size 3. More...
 
static double Norm (const Vector &a)
 Calculates the norm of vector "a". More...
 
static double StableNorm (const Vector &a)
 Calculates the norm of vector "a" while avoiding underflow and overflow. More...
 
template<class T >
static T CrossProduct (const T &a, const T &b)
 Performs the vector product of the two input vectors a,b. More...
 
template<class T1 , class T2 >
static std::enable_if< std::is_same< T1, T2 >::value, bool >::type CheckIsAlias (T1 &value1, T2 &value2)
 Checks there is aliasing. More...
 
template<class T1 , class T2 >
static std::enable_if<!std::is_same< T1, T2 >::value, bool >::type CheckIsAlias (T1 &value1, T2 &value2)
 Checks there is aliasing. More...
 
template<class T1 , class T2 , class T3 >
static void CrossProduct (T1 &c, const T2 &a, const T3 &b)
 Performs the cross product of the two input vectors a,b. More...
 
template<class T1 , class T2 , class T3 >
static void UnitCrossProduct (T1 &c, const T2 &a, const T3 &b)
 Performs the unitary cross product of the two input vectors a,b. More...
 
template<class T1 , class T2 , class T3 >
static void OrthonormalBasis (const T1 &c, T2 &a, T3 &b, const IndexType Type=0)
 This computes a orthonormal basis from a given vector (Frisvad method) More...
 
template<class T1 , class T2 , class T3 >
static void OrthonormalBasisHughesMoeller (const T1 &c, T2 &a, T3 &b)
 This computes a orthonormal basis from a given vector (Hughes Moeller method) More...
 
template<class T1 , class T2 , class T3 >
static void OrthonormalBasisFrisvad (const T1 &c, T2 &a, T3 &b)
 This computes a orthonormal basis from a given vector (Frisvad method) More...
 
template<class T1 , class T2 , class T3 >
static void OrthonormalBasisNaive (const T1 &c, T2 &a, T3 &b)
 This computes a orthonormal basis from a given vector (Naive method) More...
 
template<class T1 , class T2 >
static double VectorsAngle (const T1 &rV1, const T2 &rV2)
 Computes the angle between two vectors in 3D. More...
 
static MatrixType TensorProduct3 (const Vector &a, const Vector &b)
 Returns a matrix : A = a.tensorproduct.b. More...
 
template<class TMatrixType1 , class TMatrixType2 >
static void AddMatrix (TMatrixType1 &rDestination, const TMatrixType2 &rInputMatrix, const IndexType InitialRow, const IndexType InitialCol)
 "rInputMatrix" is ADDED to "Destination" matrix starting from InitialRow and InitialCol of the destination matrix More...
 
template<class TVectorType1 , class TVectorType2 >
static void AddVector (TVectorType1 &rDestination, const TVectorType2 &rInputVector, const IndexType InitialIndex)
 "rInputVector" is ADDED to "Destination" vector starting from InitialIndex of the destination matrix More...
 
static void SubtractMatrix (MatrixType &rDestination, const MatrixType &rInputMatrix, const IndexType InitialRow, const IndexType InitialCol)
 "rInputMatrix" is SUBTRACTED to "rDestination" matrix starting from InitialRow and InitialCol of the destination matrix More...
 
static void WriteMatrix (MatrixType &rDestination, const MatrixType &rInputMatrix, const IndexType InitialRow, const IndexType InitialCol)
 "rInputMatrix" is WRITTEN on "Destination" matrix starting from InitialRow and InitialCol of the destination matrix More...
 
static void ExpandReducedMatrix (MatrixType &rDestination, const MatrixType &rReducedMatrix, const SizeType Dimension)
 Performs the Kroneker product of the Reduced Matrix with the identity matrix of size "dimension". More...
 
static void ExpandAndAddReducedMatrix (MatrixType &rDestination, const MatrixType &rReducedMatrix, const SizeType Dimension)
 Performs the Kroneker product of the Reduced Matrix with the identity matrix of size "dimension" ADDING to the destination matrix. More...
 
static void VecAdd (Vector &rX, const double coeff, Vector &rY)
 Performs rX += coeff*rY. no check on bounds is performed. More...
 
template<class TVector , class TMatrixType = MatrixType>
static TMatrixType StressVectorToTensor (const TVector &rStressVector)
 Transforms a stess vector into a matrix. Stresses are assumed to be stored in the following way: \( [ s11, s22, s33, s12, s23, s13 ] \) for 3D case and \( [ s11, s22, s33, s12 ] \) for 2D case. \( [ s11, s22, s12 ] \) for 2D case. More...
 
template<class TVector , class TMatrixType = MatrixType>
static TMatrixType VectorToSymmetricTensor (const TVector &rVector)
 Transforms a vector into a symmetric matrix. More...
 
static int Sign (const double &ThisDataType)
 Sign function. More...
 
template<class TVector , class TMatrixType = MatrixType>
static TMatrixType StrainVectorToTensor (const TVector &rStrainVector)
 Transforms a strain vector into a matrix. Strains are assumed to be stored in the following way: \( [ e11, e22, e33, 2*e12, 2*e23, 2*e13 ] \) for 3D case and \( [ e11, e22, e33, 2*e12 ] \) for 2D case. \( [ e11, e22, 2*e12 ] \) for 2D case. More...
 
template<class TMatrixType , class TVector = Vector>
static Vector StrainTensorToVector (const TMatrixType &rStrainTensor, SizeType rSize=0)
 Transforms a given symmetric Strain Tensor to Voigt Notation: More...
 
template<class TMatrixType , class TVector = Vector>
static TVector StressTensorToVector (const TMatrixType &rStressTensor, SizeType rSize=0)
 Transforms a given symmetric Stress Tensor to Voigt Notation: More...
 
template<class TMatrixType , class TVector = Vector>
static TVector SymmetricTensorToVector (const TMatrixType &rTensor, SizeType rSize=0)
 Transforms a given symmetric Tensor to Voigt Notation: More...
 
template<class TMatrixType1 , class TMatrixType2 , class TMatrixType3 >
static void BtDBProductOperation (TMatrixType1 &rA, const TMatrixType2 &rD, const TMatrixType3 &rB)
 Calculates the product operation B'DB. More...
 
template<class TMatrixType1 , class TMatrixType2 , class TMatrixType3 >
static void BDBtProductOperation (TMatrixType1 &rA, const TMatrixType2 &rD, const TMatrixType3 &rB)
 Calculates the product operation BDB'. More...
 
template<class TMatrixType1 , class TMatrixType2 >
static bool GaussSeidelEigenSystem (const TMatrixType1 &rA, TMatrixType2 &rEigenVectorsMatrix, TMatrixType2 &rEigenValuesMatrix, const double Tolerance=1.0e-18, const SizeType MaxIterations=20)
 Calculates the eigenvectors and eigenvalues of given symmetric matrix. More...
 
template<class TIntegerType >
static TIntegerType Factorial (const TIntegerType Number)
 Calculates the Factorial of a number k, Factorial = k! More...
 
template<class TMatrixType >
static void CalculateExponentialOfMatrix (const TMatrixType &rMatrix, TMatrixType &rExponentialMatrix, const double Tolerance=1000.0 *ZeroTolerance, const SizeType MaxTerms=200)
 Calculates the exponential of a matrix. More...
 
static double DegreesToRadians (double AngleInDegrees)
 
template<SizeType TDim>
 KRATOS_DEPRECATED_MESSAGE ("Please use GaussSeidelEigenSystem() instead. Note the resulting EigenVectors matrix is transposed respect GaussSeidelEigenSystem()") static inline bool EigenSystem(const BoundedMatrix< double
 Calculates the eigenvectors and eigenvalues of given symmetric TDimxTDim matrix. More...
 
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > const double const SizeType class TMatrixType2 static inline bool MatrixSquareRoot (const TMatrixType1 &rA, TMatrixType2 &rMatrixSquareRoot, const double Tolerance=1.0e-18, const SizeType MaxIterations=20)
 

Detailed Description

template<class TDataType = double>
class Kratos::MathUtils< TDataType >

Various mathematical utilitiy functions.

Various mathematical utilitiy functions. Defines several utility functions.

Author
Riccardo Rossi
Pooyan Dadvand

Member Typedef Documentation

◆ IndexType

template<class TDataType = double>
using Kratos::MathUtils< TDataType >::IndexType = std::size_t

The index type.

◆ IndirectArrayType

template<class TDataType = double>
using Kratos::MathUtils< TDataType >::IndirectArrayType = boost::numeric::ublas::indirect_array<DenseVector<std::size_t> >

The indirect array type.

◆ MatrixType

template<class TDataType = double>
using Kratos::MathUtils< TDataType >::MatrixType = Matrix

The matrix type.

◆ SizeType

template<class TDataType = double>
using Kratos::MathUtils< TDataType >::SizeType = std::size_t

The size type.

◆ VectorType

template<class TDataType = double>
using Kratos::MathUtils< TDataType >::VectorType = Vector

The vector type.

Member Function Documentation

◆ AddMatrix()

template<class TDataType = double>
template<class TMatrixType1 , class TMatrixType2 >
static void Kratos::MathUtils< TDataType >::AddMatrix ( TMatrixType1 &  rDestination,
const TMatrixType2 &  rInputMatrix,
const IndexType  InitialRow,
const IndexType  InitialCol 
)
inlinestatic

"rInputMatrix" is ADDED to "Destination" matrix starting from InitialRow and InitialCol of the destination matrix

"Destination" is assumed to be able to contain the "input matrix" (no check is performed on the bounds)

Parameters
rDestinationThe matrix destination
rInputMatrixThe input matrix to be added
InitialRowThe initial row
InitialColThe initial column

◆ AddVector()

template<class TDataType = double>
template<class TVectorType1 , class TVectorType2 >
static void Kratos::MathUtils< TDataType >::AddVector ( TVectorType1 &  rDestination,
const TVectorType2 &  rInputVector,
const IndexType  InitialIndex 
)
inlinestatic

"rInputVector" is ADDED to "Destination" vector starting from InitialIndex of the destination matrix

"Destination" is assumed to be able to contain the "input vector" (no check is performed on the bounds)

Parameters
rDestinationThe vector destination
rInputVectorThe input vector to be added
InitialIndexThe initial index

◆ BDBtProductOperation()

template<class TDataType = double>
template<class TMatrixType1 , class TMatrixType2 , class TMatrixType3 >
static void Kratos::MathUtils< TDataType >::BDBtProductOperation ( TMatrixType1 &  rA,
const TMatrixType2 &  rD,
const TMatrixType3 &  rB 
)
inlinestatic

Calculates the product operation BDB'.

Parameters
rAThe resulting matrix
rDThe "center" matrix
rBThe matrices to be transposed
Template Parameters
TMatrixType1The type of matrix considered (1)
TMatrixType2The type of matrix considered (2)
TMatrixType3The type of matrix considered (3)

◆ BtDBProductOperation()

template<class TDataType = double>
template<class TMatrixType1 , class TMatrixType2 , class TMatrixType3 >
static void Kratos::MathUtils< TDataType >::BtDBProductOperation ( TMatrixType1 &  rA,
const TMatrixType2 &  rD,
const TMatrixType3 &  rB 
)
inlinestatic

Calculates the product operation B'DB.

Parameters
rAThe resulting matrix
rDThe "center" matrix
rBThe matrices to be transposed
Template Parameters
TMatrixType1The type of matrix considered (1)
TMatrixType2The type of matrix considered (2)
TMatrixType3The type of matrix considered (3)

◆ CalculateExponentialOfMatrix()

template<class TDataType = double>
template<class TMatrixType >
static void Kratos::MathUtils< TDataType >::CalculateExponentialOfMatrix ( const TMatrixType &  rMatrix,
TMatrixType &  rExponentialMatrix,
const double  Tolerance = 1000.0*ZeroTolerance,
const SizeType  MaxTerms = 200 
)
inlinestatic

Calculates the exponential of a matrix.

see https://mathworld.wolfram.com/MatrixExponential.html

Template Parameters
rMatrixthe matrix A of which exp is calculated
rExponentialMatrixexp(A)

◆ CheckConditionNumber()

template<class TDataType = double>
template<class TMatrix1 , class TMatrix2 >
static bool Kratos::MathUtils< TDataType >::CheckConditionNumber ( const TMatrix1 &  rInputMatrix,
TMatrix2 &  rInvertedMatrix,
const double  Tolerance = std::numeric_limits<double>::epsilon(),
const bool  ThrowError = true 
)
inlinestatic

This method checks the condition number of amtrix.

Parameters
rInputMatrixIs the input matrix (unchanged at output)
rInvertedMatrixIs the inverse of the input matrix
ToleranceThe maximum tolerance considered

◆ CheckIsAlias() [1/2]

template<class TDataType = double>
template<class T1 , class T2 >
static std::enable_if<std::is_same<T1, T2>::value, bool>::type Kratos::MathUtils< TDataType >::CheckIsAlias ( T1 &  value1,
T2 &  value2 
)
inlinestatic

Checks there is aliasing.

Parameters
value1The first value
value2The second value

◆ CheckIsAlias() [2/2]

template<class TDataType = double>
template<class T1 , class T2 >
static std::enable_if<!std::is_same<T1, T2>::value, bool>::type Kratos::MathUtils< TDataType >::CheckIsAlias ( T1 &  value1,
T2 &  value2 
)
inlinestatic

Checks there is aliasing.

Parameters
value1The first value
value2The second value

◆ Cofactor()

template<class TDataType = double>
template<class TMatrixType >
static double Kratos::MathUtils< TDataType >::Cofactor ( const TMatrixType &  rMat,
IndexType  i,
IndexType  j 
)
inlinestatic

Calculates the cofactor.

Parameters
rMatThe matrix to calculate
iThe index i
jThe index j
Returns
The cofactor of the matrix

◆ CofactorMatrix()

template<class TDataType = double>
template<class TMatrixType >
static MatrixType Kratos::MathUtils< TDataType >::CofactorMatrix ( const TMatrixType &  rMat)
inlinestatic

Calculates the cofactor matrix.

Parameters
rMatThe matrix to calculate
Returns
The cofactor matrix

◆ CrossProduct() [1/2]

template<class TDataType = double>
template<class T >
static T Kratos::MathUtils< TDataType >::CrossProduct ( const T &  a,
const T &  b 
)
inlinestatic

Performs the vector product of the two input vectors a,b.

a,b are assumed to be of size 3 (no check is performed on vector sizes)

Parameters
aFirst input vector
bSecond input vector
Returns
The resulting vector

◆ CrossProduct() [2/2]

template<class TDataType = double>
template<class T1 , class T2 , class T3 >
static void Kratos::MathUtils< TDataType >::CrossProduct ( T1 &  c,
const T2 &  a,
const T3 &  b 
)
inlinestatic

Performs the cross product of the two input vectors a,b.

a,b are assumed to be of size 3 (check is only performed on vector sizes in debug mode)

Parameters
aFirst input vector
bSecond input vector
cThe resulting vector

◆ DegreesToRadians()

template<class TDataType = double>
static double Kratos::MathUtils< TDataType >::DegreesToRadians ( double  AngleInDegrees)
inlinestatic

◆ Det()

template<class TDataType = double>
template<class TMatrixType >
static double Kratos::MathUtils< TDataType >::Det ( const TMatrixType &  rA)
inlinestatic

Calculates the determinant of a matrix of a square matrix of any size (no check performed on release mode)

Parameters
rAIs the input matrix
Returns
The determinant of any size matrix

◆ Det2()

template<class TDataType = double>
template<class TMatrixType >
static double Kratos::MathUtils< TDataType >::Det2 ( const TMatrixType &  rA)
inlinestatic

Calculates the determinant of a matrix of dimension 2x2 (no check performed on matrix size)

Parameters
rAIs the input matrix
Returns
The determinant of the 2x2 matrix

◆ Det3()

template<class TDataType = double>
template<class TMatrixType >
static double Kratos::MathUtils< TDataType >::Det3 ( const TMatrixType &  rA)
inlinestatic

Calculates the determinant of a matrix of dimension 3*3 (no check performed on matrix size)

Parameters
rAIs the input matrix
Returns
The determinant of the 3x3 matrix

◆ Det4()

template<class TDataType = double>
template<class TMatrixType >
static double Kratos::MathUtils< TDataType >::Det4 ( const TMatrixType &  rA)
inlinestatic

Calculates the determinant of a matrix of dimension 4*4 (no check performed on matrix size)

Parameters
rAIs the input matrix
Returns
The determinant of the 4x4 matrix

◆ Dot()

template<class TDataType = double>
static double Kratos::MathUtils< TDataType >::Dot ( const Vector rFirstVector,
const Vector rSecondVector 
)
inlinestatic

Performs the dot product of two vectors of arbitrary size.

No check performed on vector sizes

Parameters
rFirstVectorFirst input vector
rSecondVectorSecond input vector
Returns
The resulting norm

◆ Dot3()

template<class TDataType = double>
static double Kratos::MathUtils< TDataType >::Dot3 ( const Vector a,
const Vector b 
)
inlinestatic

Performs the dot product of two vectors of dimension 3.

No check performed on vector sizes

Parameters
aFirst input vector
bSecond input vector
Returns
The resulting norm

◆ ExpandAndAddReducedMatrix()

template<class TDataType = double>
static void Kratos::MathUtils< TDataType >::ExpandAndAddReducedMatrix ( MatrixType rDestination,
const MatrixType rReducedMatrix,
const SizeType  Dimension 
)
inlinestatic

Performs the Kroneker product of the Reduced Matrix with the identity matrix of size "dimension" ADDING to the destination matrix.

Parameters
rDestinationThe matric destination
rReducedMatrixThe reduced matrix to be added
DimensionThe dimension where we work

◆ ExpandReducedMatrix()

template<class TDataType = double>
static void Kratos::MathUtils< TDataType >::ExpandReducedMatrix ( MatrixType rDestination,
const MatrixType rReducedMatrix,
const SizeType  Dimension 
)
inlinestatic

Performs the Kroneker product of the Reduced Matrix with the identity matrix of size "dimension".

Parameters
rDestinationThe matric destination
rReducedMatrixThe reduced matrix to be computed
DimensionThe dimension where we work

◆ Factorial()

template<class TDataType = double>
template<class TIntegerType >
static TIntegerType Kratos::MathUtils< TDataType >::Factorial ( const TIntegerType  Number)
inlinestatic

Calculates the Factorial of a number k, Factorial = k!

Template Parameters
NumberThe number of which the Factorial is computed

◆ GaussSeidelEigenSystem()

template<class TDataType = double>
template<class TMatrixType1 , class TMatrixType2 >
static bool Kratos::MathUtils< TDataType >::GaussSeidelEigenSystem ( const TMatrixType1 &  rA,
TMatrixType2 &  rEigenVectorsMatrix,
TMatrixType2 &  rEigenValuesMatrix,
const double  Tolerance = 1.0e-18,
const SizeType  MaxIterations = 20 
)
inlinestatic

Calculates the eigenvectors and eigenvalues of given symmetric matrix.

The eigenvectors and eigenvalues are calculated using the iterative Gauss-Seidel-method. The resulting decomposition is LDL'

Note
See https://en.wikipedia.org/wiki/Gauss%E2%80%93Seidel_method
Parameters
rAThe given symmetric matrix the eigenvectors are to be calculated.
rEigenVectorsMatrixThe result matrix (will be overwritten with the eigenvectors)
rEigenValuesMatrixThe result diagonal matrix with the eigenvalues
ToleranceThe largest value considered to be zero
MaxIterationsMaximum number of iterations
Template Parameters
TMatrixType1The type of matrix considered (1)
TMatrixType2The type of matrix considered (2)

◆ GeneralizedDet()

template<class TDataType = double>
template<class TMatrixType >
static double Kratos::MathUtils< TDataType >::GeneralizedDet ( const TMatrixType &  rA)
inlinestatic

Calculates the determinant of any matrix (no check performed on matrix size)

Parameters
rAIs the input matrix
Returns
The determinant of the 2x2 matrix

◆ GeneralizedInvertMatrix()

template<class TDataType = double>
template<class TMatrix1 , class TMatrix2 >
static void Kratos::MathUtils< TDataType >::GeneralizedInvertMatrix ( const TMatrix1 &  rInputMatrix,
TMatrix2 &  rInvertedMatrix,
double rInputMatrixDet,
const double  Tolerance = ZeroTolerance 
)
inlinestatic

It inverts non square matrices (https://en.wikipedia.org/wiki/Inverse_element#Matrices)

Parameters
rInputMatrixIs the input matrix (unchanged at output)
rInvertedMatrixIs the inverse of the input matrix
rInputMatrixDetIs the determinant of the input matrix
ToleranceThe maximum tolerance considered
Template Parameters
TMatrix1The type of the input matrix
TMatrix2Is the type of the output matrix

◆ GetZeroTolerance()

template<class TDataType = double>
static double Kratos::MathUtils< TDataType >::GetZeroTolerance ( )
inlinestatic

This function returns the machine precision.

Returns
The corresponding epsilon for the double

◆ Heron()

template<class TDataType = double>
template<bool TCheck>
static double Kratos::MathUtils< TDataType >::Heron ( double  a,
double  b,
double  c 
)
inlinestatic

In geometry, Heron's formula (sometimes called Hero's formula), named after Hero of Alexandria, gives the area of a triangle by requiring no arbitrary choice of side as base or vertex as origin, contrary to other formulas for the area of a triangle, such as half the base times the height or half the norm of a cross product of two sides.

Parameters
aFirst length
bSecond length
cThird length
Returns
Heron solution: Heron's formula states that the area of a triangle whose sides have lengths a, b, and c

◆ InvertMatrix() [1/2]

template<class TDataType = double>
template<SizeType TDim>
static BoundedMatrix<double, TDim, TDim> Kratos::MathUtils< TDataType >::InvertMatrix ( const BoundedMatrix< double, TDim, TDim > &  rInputMatrix,
double rInputMatrixDet,
const double  Tolerance = ZeroTolerance 
)
inlinestatic

Calculates the inverse of a 2x2, 3x3 and 4x4 matrices (using bounded matrix for performance)

Parameters
rInputMatrixThe matrix to invert
rInputMatrixDetThe determinant of the matrix
ToleranceThe maximum tolerance considered
Returns
InvertMatrix: The inverted matrix

◆ InvertMatrix() [2/2]

template<class TDataType = double>
template<class TMatrix1 , class TMatrix2 >
static void Kratos::MathUtils< TDataType >::InvertMatrix ( const TMatrix1 &  rInputMatrix,
TMatrix2 &  rInvertedMatrix,
double rInputMatrixDet,
const double  Tolerance = ZeroTolerance 
)
inlinestatic

It inverts matrices of order 2, 3 and 4.

Parameters
rInputMatrixIs the input matrix (unchanged at output)
rInvertedMatrixIs the inverse of the input matrix
rInputMatrixDetIs the determinant of the input matrix
Template Parameters
TMatrix1The type of the input matrix
TMatrix2Is the type of the output matrix

◆ InvertMatrix2()

template<class TDataType = double>
template<class TMatrix1 , class TMatrix2 >
static void Kratos::MathUtils< TDataType >::InvertMatrix2 ( const TMatrix1 &  rInputMatrix,
TMatrix2 &  rInvertedMatrix,
double rInputMatrixDet 
)
inlinestatic

It inverts matrices of order 2.

Parameters
rInputMatrixIs the input matrix (unchanged at output)
rInvertedMatrixIs the inverse of the input matrix
rInputMatrixDetIs the determinant of the input matrix

◆ InvertMatrix3()

template<class TDataType = double>
template<class TMatrix1 , class TMatrix2 >
static void Kratos::MathUtils< TDataType >::InvertMatrix3 ( const TMatrix1 &  rInputMatrix,
TMatrix2 &  rInvertedMatrix,
double rInputMatrixDet 
)
inlinestatic

It inverts matrices of order 3.

Parameters
rInputMatrixIs the input matrix (unchanged at output)
rInvertedMatrixIs the inverse of the input matrix
rInputMatrixDetIs the determinant of the input matrix

◆ InvertMatrix4()

template<class TDataType = double>
template<class TMatrix1 , class TMatrix2 >
static void Kratos::MathUtils< TDataType >::InvertMatrix4 ( const TMatrix1 &  rInputMatrix,
TMatrix2 &  rInvertedMatrix,
double rInputMatrixDet 
)
inlinestatic

It inverts matrices of order 4.

Parameters
rInputMatrixIs the input matrix (unchanged at output)
rInvertedMatrixIs the inverse of the input matrix
rInputMatrixDetIs the determinant of the input matrix

◆ KRATOS_DEPRECATED_MESSAGE()

template<class TDataType = double>
template<SizeType TDim>
Kratos::MathUtils< TDataType >::KRATOS_DEPRECATED_MESSAGE ( "Please use GaussSeidelEigenSystem() instead. Note the resulting EigenVectors matrix is transposed respect GaussSeidelEigenSystem()"  ) const

Calculates the eigenvectors and eigenvalues of given symmetric TDimxTDim matrix.

The eigenvectors and eigenvalues are calculated using the iterative Gauss-Seidel-method. The resulting decomposition is L'DL

Note
See https://en.wikipedia.org/wiki/Gauss%E2%80%93Seidel_method
Parameters
AThe given symmetric matrix the eigenvectors are to be calculated.
rEigenVectorsMatrixThe result matrix (will be overwritten with the eigenvectors)
rEigenValuesMatrixThe result diagonal matrix with the eigenvalues
ToleranceThe largest value considered to be zero
MaxIterationsMaximum number of iterations
Template Parameters
TDimThe size of the bounded matrix
Warning
This method is deprecated. Will be removed soon

◆ MatrixSquareRoot()

template<class TDataType = double>
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim> const double const SizeType class TMatrixType2 static inline bool Kratos::MathUtils< TDataType >::MatrixSquareRoot ( const TMatrixType1 &  rA,
TMatrixType2 &  rMatrixSquareRoot,
const double  Tolerance = 1.0e-18,
const SizeType  MaxIterations = 20 
)
inline

◆ Norm()

template<class TDataType = double>
static double Kratos::MathUtils< TDataType >::Norm ( const Vector a)
inlinestatic

Calculates the norm of vector "a".

Parameters
aInput vector
Returns
The resulting norm

◆ Norm3()

template<class TDataType = double>
template<class TVectorType >
static double Kratos::MathUtils< TDataType >::Norm3 ( const TVectorType &  a)
inlinestatic

Calculates the norm of vector "a" which is assumed to be of size 3.

No check is performed on the vector's size

Parameters
aInput vector
Returns
The resulting norm

◆ OrthonormalBasis()

template<class TDataType = double>
template<class T1 , class T2 , class T3 >
static void Kratos::MathUtils< TDataType >::OrthonormalBasis ( const T1 &  c,
T2 &  a,
T3 &  b,
const IndexType  Type = 0 
)
inlinestatic

This computes a orthonormal basis from a given vector (Frisvad method)

Parameters
cThe input vector
aFirst resulting vector
bSecond resulting vector
TypeThe type of method employed, 0 is HughesMoeller, 1 is Frisvad and otherwise Naive

◆ OrthonormalBasisFrisvad()

template<class TDataType = double>
template<class T1 , class T2 , class T3 >
static void Kratos::MathUtils< TDataType >::OrthonormalBasisFrisvad ( const T1 &  c,
T2 &  a,
T3 &  b 
)
inlinestatic

This computes a orthonormal basis from a given vector (Frisvad method)

Parameters
cThe input vector
aFirst resulting vector
bSecond resulting vector
Note
Orthonormal basis taken from: http://orbit.dtu.dk/files/126824972/onb_frisvad_jgt2012_v2.pdf

◆ OrthonormalBasisHughesMoeller()

template<class TDataType = double>
template<class T1 , class T2 , class T3 >
static void Kratos::MathUtils< TDataType >::OrthonormalBasisHughesMoeller ( const T1 &  c,
T2 &  a,
T3 &  b 
)
inlinestatic

This computes a orthonormal basis from a given vector (Hughes Moeller method)

Parameters
cThe input vector
aFirst resulting vector
bSecond resulting vector
Note
Orthonormal basis taken from: http://orbit.dtu.dk/files/126824972/onb_frisvad_jgt2012_v2.pdf

◆ OrthonormalBasisNaive()

template<class TDataType = double>
template<class T1 , class T2 , class T3 >
static void Kratos::MathUtils< TDataType >::OrthonormalBasisNaive ( const T1 &  c,
T2 &  a,
T3 &  b 
)
inlinestatic

This computes a orthonormal basis from a given vector (Naive method)

Parameters
cThe input vector
aFirst resulting vector
bSecond resulting vector
Note
Orthonormal basis taken from: http://orbit.dtu.dk/files/126824972/onb_frisvad_jgt2012_v2.pdf

◆ Sign()

template<class TDataType = double>
static int Kratos::MathUtils< TDataType >::Sign ( const double ThisDataType)
inlinestatic

Sign function.

Parameters
ThisDataTypeThe value to extract the sign
Returns
The sign of the value

◆ Solve()

template<class TDataType >
void Kratos::MathUtils< TDataType >::Solve ( MatrixType  A,
VectorType rX,
const VectorType rB 
)
static

This function is designed to be called when a dense linear system is needed to be solved.

Parameters
ASystem matrix
rXSolution vector. it's also the initial guess for iterative linear solvers.
rBRight hand side vector.

◆ StableNorm()

template<class TDataType = double>
static double Kratos::MathUtils< TDataType >::StableNorm ( const Vector a)
inlinestatic

Calculates the norm of vector "a" while avoiding underflow and overflow.

Parameters
aInput vector
Returns
The resulting norm
See also
http://www.netlib.org/lapack/explore-html/da/d7f/dnrm2_8f_source.html

◆ StrainTensorToVector()

template<class TDataType = double>
template<class TMatrixType , class TVector = Vector>
static Vector Kratos::MathUtils< TDataType >::StrainTensorToVector ( const TMatrixType &  rStrainTensor,
SizeType  rSize = 0 
)
inlinestatic

Transforms a given symmetric Strain Tensor to Voigt Notation:

The following cases:

  • In the 3D case: from a second order tensor (3*3) Matrix to a corresponing (6*1) Vector \( [ e11, e22, e33, 2*e12, 2*e23, 2*e13 ] \) for 3D case and
  • In the 2D case: from a second order tensor (3*3) Matrix to a corresponing (4*1) Vector \( [ e11, e22, e33, 2*e12 ] \) fir 2D case.
  • In the 2D case: from a second order tensor (2*2) Matrix to a corresponing (3*1) Vector \( [ e11, e22, 2*e12 ] \) fir 2D case.
    Parameters
    rStrainTensorthe given symmetric second order strain tensor
    Returns
    the corresponding strain tensor in vector form
    Template Parameters
    TMatrixTypeThe matrix type considered
    TVectorThe vector returning type

◆ StrainVectorToTensor()

template<class TDataType = double>
template<class TVector , class TMatrixType = MatrixType>
static TMatrixType Kratos::MathUtils< TDataType >::StrainVectorToTensor ( const TVector &  rStrainVector)
inlinestatic

Transforms a strain vector into a matrix. Strains are assumed to be stored in the following way: \( [ e11, e22, e33, 2*e12, 2*e23, 2*e13 ] \) for 3D case and \( [ e11, e22, e33, 2*e12 ] \) for 2D case. \( [ e11, e22, 2*e12 ] \) for 2D case.

Hence the deviatoric components of the strain vector are divided by 2 while they are stored into the matrix

Parameters
rStrainVectorthe given strain vector
Returns
the corresponding strain tensor in matrix form
Template Parameters
TVectorThe vector type considered
TMatrixTypeThe matrix returning type

◆ StressTensorToVector()

template<class TDataType = double>
template<class TMatrixType , class TVector = Vector>
static TVector Kratos::MathUtils< TDataType >::StressTensorToVector ( const TMatrixType &  rStressTensor,
SizeType  rSize = 0 
)
inlinestatic

Transforms a given symmetric Stress Tensor to Voigt Notation:

Components are assumed to be stored in the following way: \( [ s11, s22, s33, s12, s23, s13 ] \) for 3D case and \( [ s11, s22, s33, s12 ] \) for 2D case. \( [ s11, s22, s12 ] \) for 2D case. In the 3D case: from a second order tensor (3*3) Matrix to a corresponing (6*1) Vector In the 3D case: from a second order tensor (3*3) Matrix to a corresponing (4*1) Vector In the 2D case: from a second order tensor (2*2) Matrix to a corresponing (3*1) Vector

Parameters
rStressTensorthe given symmetric second order stress tensor
Returns
the corresponding stress tensor in vector form
Template Parameters
TMatrixTypeThe matrix type considered
TVectorThe vector returning type

◆ StressVectorToTensor()

template<class TDataType = double>
template<class TVector , class TMatrixType = MatrixType>
static TMatrixType Kratos::MathUtils< TDataType >::StressVectorToTensor ( const TVector &  rStressVector)
inlinestatic

Transforms a stess vector into a matrix. Stresses are assumed to be stored in the following way: \( [ s11, s22, s33, s12, s23, s13 ] \) for 3D case and \( [ s11, s22, s33, s12 ] \) for 2D case. \( [ s11, s22, s12 ] \) for 2D case.

Parameters
rStressVectorthe given stress vector
Returns
the corresponding stress tensor in matrix form
Template Parameters
TVectorThe vector type considered
TMatrixTypeThe matrix returning type

◆ SubtractMatrix()

template<class TDataType = double>
static void Kratos::MathUtils< TDataType >::SubtractMatrix ( MatrixType rDestination,
const MatrixType rInputMatrix,
const IndexType  InitialRow,
const IndexType  InitialCol 
)
inlinestatic

"rInputMatrix" is SUBTRACTED to "rDestination" matrix starting from InitialRow and InitialCol of the destination matrix

"rDestination" is assumed to be able to contain the "input matrix" (no check is performed on the bounds)

Parameters
rDestinationThe matric destination
rInputMatrixThe input matrix to be computed
InitialRowThe initial row to compute
InitialColThe initial column to compute

◆ SymmetricTensorToVector()

template<class TDataType = double>
template<class TMatrixType , class TVector = Vector>
static TVector Kratos::MathUtils< TDataType >::SymmetricTensorToVector ( const TMatrixType &  rTensor,
SizeType  rSize = 0 
)
inlinestatic

Transforms a given symmetric Tensor to Voigt Notation:

The following cases:

  • In the 3D case: from a second order tensor (3*3) Matrix to a corresponing (6*1) Vector
  • In the 3D case: from a second order tensor (3*3) Matrix to a corresponing (4*1) Vector
  • In the 2D case: from a second order tensor (2*2) Matrix to a corresponing (3*1) Vector
    Parameters
    rTensorthe given symmetric second order stress tensor
    Returns
    the corresponding stress tensor in vector form
    Template Parameters
    TMatrixTypeThe matrix type considered
    TVectorThe vector returning type

◆ TensorProduct3()

template<class TDataType = double>
static MatrixType Kratos::MathUtils< TDataType >::TensorProduct3 ( const Vector a,
const Vector b 
)
inlinestatic

Returns a matrix : A = a.tensorproduct.b.

a,b are assumed to be of order 3, no check is performed on the size of the vectors

Parameters
aFirst input vector
bSecond input vector
Returns
Returns A = a.tensorproduct.b

◆ UnitCrossProduct()

template<class TDataType = double>
template<class T1 , class T2 , class T3 >
static void Kratos::MathUtils< TDataType >::UnitCrossProduct ( T1 &  c,
const T2 &  a,
const T3 &  b 
)
inlinestatic

Performs the unitary cross product of the two input vectors a,b.

a,b are assumed to be of size 3 (no check is performed on vector sizes)

Parameters
aFirst input vector
bSecond input vector
cThe resulting vector

◆ VecAdd()

template<class TDataType = double>
static void Kratos::MathUtils< TDataType >::VecAdd ( Vector rX,
const double  coeff,
Vector rY 
)
inlinestatic

Performs rX += coeff*rY. no check on bounds is performed.

Parameters
rXThe vector destination
rYThe vector to be added
coeffThe proportion to be added

◆ VectorsAngle()

template<class TDataType = double>
template<class T1 , class T2 >
static double Kratos::MathUtils< TDataType >::VectorsAngle ( const T1 &  rV1,
const T2 &  rV2 
)
inlinestatic

Computes the angle between two vectors in 3D.

Parameters
rV1First input vector
rV2Second input vector

◆ VectorToSymmetricTensor()

template<class TDataType = double>
template<class TVector , class TMatrixType = MatrixType>
static TMatrixType Kratos::MathUtils< TDataType >::VectorToSymmetricTensor ( const TVector &  rVector)
inlinestatic

Transforms a vector into a symmetric matrix.

Components are assumed to be stored in the following way: \( [ s11, s22, s33, s12, s23, s13 ] \) for 3D case and \( [ s11, s22, s33, s12 ] \) for 2D case. \( [ s11, s22, s12 ] \) for 2D case.

Parameters
rVectorthe given stress vector
Returns
The corresponding Tensor in matrix form
Template Parameters
TVectorThe vector type considered
TMatrixTypeThe matrix returning type

◆ WriteMatrix()

template<class TDataType = double>
static void Kratos::MathUtils< TDataType >::WriteMatrix ( MatrixType rDestination,
const MatrixType rInputMatrix,
const IndexType  InitialRow,
const IndexType  InitialCol 
)
inlinestatic

"rInputMatrix" is WRITTEN on "Destination" matrix starting from InitialRow and InitialCol of the destination matrix

"Destination" is assumed to be able to contain the "input matrix" (no check is performed on the bounds)

Warning
Destination is overwritten!!
Parameters
rDestinationThe matric destination
rrInputMatrixThe input matrix to be computed
InitialRowThe initial row to compute
InitialColThe initial column to compute

Member Data Documentation

◆ MaxIterations

template<class TDataType = double>
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim> const double const SizeType Kratos::MathUtils< TDataType >::MaxIterations
Initial value:
= 20
)
{
const BoundedMatrix<double, TDim, TDim> V_matrix = rEigenVectorsMatrix;
for(IndexType i = 0; i < TDim; ++i) {
for(IndexType j = 0; j < TDim; ++j) {
rEigenVectorsMatrix(i, j) = V_matrix(j, i);
}
}
return is_converged;
}
template<class TMatrixType1
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > const double Tolerance
Definition: math_utils.h:1715
std::size_t IndexType
The index type.
Definition: math_utils.h:78
TDim & rA
Definition: math_utils.h:1712
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > & rEigenValuesMatrix
Definition: math_utils.h:1714
TDim
Definition: math_utils.h:1712
TDim BoundedMatrix< double, TDim, TDim > BoundedMatrix< double, TDim, TDim > const double const SizeType MaxIterations
Definition: math_utils.h:1716
static bool GaussSeidelEigenSystem(const TMatrixType1 &rA, TMatrixType2 &rEigenVectorsMatrix, TMatrixType2 &rEigenValuesMatrix, const double Tolerance=1.0e-18, const SizeType MaxIterations=20)
Calculates the eigenvectors and eigenvalues of given symmetric matrix.
Definition: math_utils.h:1587
TDim BoundedMatrix< double, TDim, TDim > & rEigenVectorsMatrix
Definition: math_utils.h:1713
int j
Definition: quadrature.py:648
integer i
Definition: TensorModule.f:17

◆ rA

template<class TDataType = double>
TDim& Kratos::MathUtils< TDataType >::rA

◆ rEigenValuesMatrix

template<class TDataType = double>
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim>& Kratos::MathUtils< TDataType >::rEigenValuesMatrix

◆ rEigenVectorsMatrix

template<class TDataType = double>
TDim BoundedMatrix<double, TDim, TDim>& Kratos::MathUtils< TDataType >::rEigenVectorsMatrix

◆ TDim

template<class TDataType = double>
Kratos::MathUtils< TDataType >::TDim

◆ Tolerance

template<class TDataType = double>
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim> const double Kratos::MathUtils< TDataType >::Tolerance = 1.0e-18

◆ ZeroTolerance

template<class TDataType = double>
constexpr double Kratos::MathUtils< TDataType >::ZeroTolerance = std::numeric_limits<double>::epsilon()
staticconstexpr

The machine precision.


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