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.
|
Various mathematical utilitiy functions. More...
#include <math_utils.h>
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 | |
TDim & | rA |
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, TDim > | InvertMatrix (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) |
Various mathematical utilitiy functions.
Various mathematical utilitiy functions. Defines several utility functions.
using Kratos::MathUtils< TDataType >::IndexType = std::size_t |
The index type.
using Kratos::MathUtils< TDataType >::IndirectArrayType = boost::numeric::ublas::indirect_array<DenseVector<std::size_t> > |
The indirect array type.
using Kratos::MathUtils< TDataType >::MatrixType = Matrix |
The matrix type.
using Kratos::MathUtils< TDataType >::SizeType = std::size_t |
The size type.
using Kratos::MathUtils< TDataType >::VectorType = Vector |
The vector type.
|
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)
rDestination | The matrix destination |
rInputMatrix | The input matrix to be added |
InitialRow | The initial row |
InitialCol | The initial column |
|
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)
rDestination | The vector destination |
rInputVector | The input vector to be added |
InitialIndex | The initial index |
|
inlinestatic |
Calculates the product operation BDB'.
rA | The resulting matrix |
rD | The "center" matrix |
rB | The matrices to be transposed |
TMatrixType1 | The type of matrix considered (1) |
TMatrixType2 | The type of matrix considered (2) |
TMatrixType3 | The type of matrix considered (3) |
|
inlinestatic |
Calculates the product operation B'DB.
rA | The resulting matrix |
rD | The "center" matrix |
rB | The matrices to be transposed |
TMatrixType1 | The type of matrix considered (1) |
TMatrixType2 | The type of matrix considered (2) |
TMatrixType3 | The type of matrix considered (3) |
|
inlinestatic |
Calculates the exponential of a matrix.
see https://mathworld.wolfram.com/MatrixExponential.html
rMatrix | the matrix A of which exp is calculated |
rExponentialMatrix | exp(A) |
|
inlinestatic |
This method checks the condition number of amtrix.
rInputMatrix | Is the input matrix (unchanged at output) |
rInvertedMatrix | Is the inverse of the input matrix |
Tolerance | The maximum tolerance considered |
|
inlinestatic |
Checks there is aliasing.
value1 | The first value |
value2 | The second value |
|
inlinestatic |
Checks there is aliasing.
value1 | The first value |
value2 | The second value |
|
inlinestatic |
Calculates the cofactor.
rMat | The matrix to calculate |
i | The index i |
j | The index j |
|
inlinestatic |
Calculates the cofactor matrix.
rMat | The matrix to calculate |
|
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)
a | First input vector |
b | Second input vector |
|
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)
a | First input vector |
b | Second input vector |
c | The resulting vector |
|
inlinestatic |
|
inlinestatic |
Calculates the determinant of a matrix of a square matrix of any size (no check performed on release mode)
rA | Is the input matrix |
|
inlinestatic |
Calculates the determinant of a matrix of dimension 2x2 (no check performed on matrix size)
rA | Is the input matrix |
|
inlinestatic |
Calculates the determinant of a matrix of dimension 3*3 (no check performed on matrix size)
rA | Is the input matrix |
|
inlinestatic |
Calculates the determinant of a matrix of dimension 4*4 (no check performed on matrix size)
rA | Is the input matrix |
|
inlinestatic |
Performs the dot product of two vectors of arbitrary size.
No check performed on vector sizes
rFirstVector | First input vector |
rSecondVector | Second input vector |
|
inlinestatic |
Performs the dot product of two vectors of dimension 3.
No check performed on vector sizes
a | First input vector |
b | Second input vector |
|
inlinestatic |
Performs the Kroneker product of the Reduced Matrix with the identity matrix of size "dimension" ADDING to the destination matrix.
rDestination | The matric destination |
rReducedMatrix | The reduced matrix to be added |
Dimension | The dimension where we work |
|
inlinestatic |
Performs the Kroneker product of the Reduced Matrix with the identity matrix of size "dimension".
rDestination | The matric destination |
rReducedMatrix | The reduced matrix to be computed |
Dimension | The dimension where we work |
|
inlinestatic |
Calculates the Factorial of a number k, Factorial = k!
Number | The number of which the Factorial is computed |
|
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'
rA | The given symmetric matrix the eigenvectors are to be calculated. |
rEigenVectorsMatrix | The result matrix (will be overwritten with the eigenvectors) |
rEigenValuesMatrix | The result diagonal matrix with the eigenvalues |
Tolerance | The largest value considered to be zero |
MaxIterations | Maximum number of iterations |
TMatrixType1 | The type of matrix considered (1) |
TMatrixType2 | The type of matrix considered (2) |
|
inlinestatic |
Calculates the determinant of any matrix (no check performed on matrix size)
rA | Is the input matrix |
|
inlinestatic |
It inverts non square matrices (https://en.wikipedia.org/wiki/Inverse_element#Matrices)
rInputMatrix | Is the input matrix (unchanged at output) |
rInvertedMatrix | Is the inverse of the input matrix |
rInputMatrixDet | Is the determinant of the input matrix |
Tolerance | The maximum tolerance considered |
TMatrix1 | The type of the input matrix |
TMatrix2 | Is the type of the output matrix |
|
inlinestatic |
This function returns the machine precision.
|
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.
a | First length |
b | Second length |
c | Third length |
|
inlinestatic |
Calculates the inverse of a 2x2, 3x3 and 4x4 matrices (using bounded matrix for performance)
rInputMatrix | The matrix to invert |
rInputMatrixDet | The determinant of the matrix |
Tolerance | The maximum tolerance considered |
|
inlinestatic |
It inverts matrices of order 2, 3 and 4.
rInputMatrix | Is the input matrix (unchanged at output) |
rInvertedMatrix | Is the inverse of the input matrix |
rInputMatrixDet | Is the determinant of the input matrix |
TMatrix1 | The type of the input matrix |
TMatrix2 | Is the type of the output matrix |
|
inlinestatic |
It inverts matrices of order 2.
rInputMatrix | Is the input matrix (unchanged at output) |
rInvertedMatrix | Is the inverse of the input matrix |
rInputMatrixDet | Is the determinant of the input matrix |
|
inlinestatic |
It inverts matrices of order 3.
rInputMatrix | Is the input matrix (unchanged at output) |
rInvertedMatrix | Is the inverse of the input matrix |
rInputMatrixDet | Is the determinant of the input matrix |
|
inlinestatic |
It inverts matrices of order 4.
rInputMatrix | Is the input matrix (unchanged at output) |
rInvertedMatrix | Is the inverse of the input matrix |
rInputMatrixDet | Is the determinant of the input matrix |
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
A | The given symmetric matrix the eigenvectors are to be calculated. |
rEigenVectorsMatrix | The result matrix (will be overwritten with the eigenvectors) |
rEigenValuesMatrix | The result diagonal matrix with the eigenvalues |
Tolerance | The largest value considered to be zero |
MaxIterations | Maximum number of iterations |
TDim | The size of the bounded matrix |
|
inline |
|
inlinestatic |
Calculates the norm of vector "a".
a | Input vector |
|
inlinestatic |
Calculates the norm of vector "a" which is assumed to be of size 3.
No check is performed on the vector's size
a | Input vector |
|
inlinestatic |
This computes a orthonormal basis from a given vector (Frisvad method)
c | The input vector |
a | First resulting vector |
b | Second resulting vector |
Type | The type of method employed, 0 is HughesMoeller, 1 is Frisvad and otherwise Naive |
|
inlinestatic |
This computes a orthonormal basis from a given vector (Frisvad method)
c | The input vector |
a | First resulting vector |
b | Second resulting vector |
|
inlinestatic |
This computes a orthonormal basis from a given vector (Hughes Moeller method)
c | The input vector |
a | First resulting vector |
b | Second resulting vector |
|
inlinestatic |
This computes a orthonormal basis from a given vector (Naive method)
c | The input vector |
a | First resulting vector |
b | Second resulting vector |
|
inlinestatic |
Sign function.
ThisDataType | The value to extract the sign |
|
static |
This function is designed to be called when a dense linear system is needed to be solved.
A | System matrix |
rX | Solution vector. it's also the initial guess for iterative linear solvers. |
rB | Right hand side vector. |
|
inlinestatic |
Calculates the norm of vector "a" while avoiding underflow and overflow.
a | Input vector |
|
inlinestatic |
Transforms a given symmetric Strain Tensor to Voigt Notation:
The following cases:
rStrainTensor | the given symmetric second order strain tensor |
TMatrixType | The matrix type considered |
TVector | The vector returning type |
|
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
rStrainVector | the given strain vector |
TVector | The vector type considered |
TMatrixType | The matrix returning type |
|
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
rStressTensor | the given symmetric second order stress tensor |
TMatrixType | The matrix type considered |
TVector | The vector returning type |
|
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.
rStressVector | the given stress vector |
TVector | The vector type considered |
TMatrixType | The matrix returning type |
|
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)
rDestination | The matric destination |
rInputMatrix | The input matrix to be computed |
InitialRow | The initial row to compute |
InitialCol | The initial column to compute |
|
inlinestatic |
Transforms a given symmetric Tensor to Voigt Notation:
The following cases:
rTensor | the given symmetric second order stress tensor |
TMatrixType | The matrix type considered |
TVector | The vector returning type |
|
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
a | First input vector |
b | Second input vector |
|
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)
a | First input vector |
b | Second input vector |
c | The resulting vector |
|
inlinestatic |
Performs rX += coeff*rY. no check on bounds is performed.
rX | The vector destination |
rY | The vector to be added |
coeff | The proportion to be added |
|
inlinestatic |
Computes the angle between two vectors in 3D.
rV1 | First input vector |
rV2 | Second input vector |
|
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.
rVector | the given stress vector |
TVector | The vector type considered |
TMatrixType | The matrix returning type |
|
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)
rDestination | The matric destination |
rrInputMatrix | The input matrix to be computed |
InitialRow | The initial row to compute |
InitialCol | The initial column to compute |
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim> const double const SizeType Kratos::MathUtils< TDataType >::MaxIterations |
TDim& Kratos::MathUtils< TDataType >::rA |
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim>& Kratos::MathUtils< TDataType >::rEigenValuesMatrix |
TDim BoundedMatrix<double, TDim, TDim>& Kratos::MathUtils< TDataType >::rEigenVectorsMatrix |
Kratos::MathUtils< TDataType >::TDim |
TDim BoundedMatrix<double, TDim, TDim> BoundedMatrix<double, TDim, TDim> const double Kratos::MathUtils< TDataType >::Tolerance = 1.0e-18 |
|
staticconstexpr |
The machine precision.