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.
Static Public Member Functions | List of all members
Kratos::GeoMechanicsMathUtilities< TDataType > Class Template Reference

#include <math_utilities.hpp>

Collaboration diagram for Kratos::GeoMechanicsMathUtilities< TDataType >:

Public Types

type definitions
typedef Matrix MatrixType
 
typedef Vector VectorType
 
typedef unsigned int IndexType
 
typedef unsigned int SizeType
 
typedef MathUtils< TDataType > MathUtilsType
 
typedef DenseVector< VectorSecond_Order_Tensor
 
typedef DenseVector< Second_Order_TensorThird_Order_Tensor
 
typedef DenseVector< DenseVector< Matrix > > Fourth_Order_Tensor
 
typedef matrix< Second_Order_TensorMatrix_Second_Tensor
 

Static Public Member Functions

static bool CardanoFormula (double a, double b, double c, double d, Vector &solution)
 
static Vector EigenValues (const Matrix &A, double tolerance, double zero)
 
static Vector EigenValuesDirectMethod (const Matrix &A)
 
static void QRFactorization (const MatrixType &A, MatrixType &Q, MatrixType &R)
 
static void EigenVectors (const MatrixType &A, MatrixType &vectors, VectorType &lambda, double zero_tolerance=1e-9, int max_iterations=10)
 
static MatrixType StrainVectorToTensor (const VectorType &Strains)
 
static Vector TensorToStrainVector (const Matrix &Tensor)
 
static double NormTensor (Matrix &Tensor)
 
static void VectorToTensor (const Vector &Stress, Matrix &Tensor)
 
static void TensorToVector (const Matrix &Tensor, Vector &Vector)
 
static void TensorToMatrix (Fourth_Order_Tensor &Tensor, Matrix &Matrix)
 
static void MatrixToTensor (MatrixType &A, std::vector< std::vector< Matrix > > &Tensor)
 
static void MatrixToTensor (MatrixType &A, array_1d< double, 81 > &Tensor)
 
static void TensorToMatrix (std::vector< std::vector< Matrix > > &Tensor, Matrix &Matrix)
 
static void TensorToMatrix (const array_1d< double, 81 > &Tensor, Matrix &Matrix)
 
static void DeviatoricUnity (std::vector< std::vector< Matrix > > &Unity)
 
static void DeviatoricUnity (array_1d< double, 81 > &Unity)
 
static bool Clipping (std::vector< Point * > &clipping_points, std::vector< Point * > &subjected_points, std::vector< Point * > &result_points, double tolerance)
 

Member Typedef Documentation

◆ Fourth_Order_Tensor

template<class TDataType >
typedef DenseVector<DenseVector<Matrix> > Kratos::GeoMechanicsMathUtilities< TDataType >::Fourth_Order_Tensor

◆ IndexType

template<class TDataType >
typedef unsigned int Kratos::GeoMechanicsMathUtilities< TDataType >::IndexType

◆ MathUtilsType

template<class TDataType >
typedef MathUtils<TDataType> Kratos::GeoMechanicsMathUtilities< TDataType >::MathUtilsType

◆ Matrix_Second_Tensor

template<class TDataType >
typedef matrix<Second_Order_Tensor> Kratos::GeoMechanicsMathUtilities< TDataType >::Matrix_Second_Tensor

◆ MatrixType

template<class TDataType >
typedef Matrix Kratos::GeoMechanicsMathUtilities< TDataType >::MatrixType

◆ Second_Order_Tensor

template<class TDataType >
typedef DenseVector<Vector> Kratos::GeoMechanicsMathUtilities< TDataType >::Second_Order_Tensor

◆ SizeType

template<class TDataType >
typedef unsigned int Kratos::GeoMechanicsMathUtilities< TDataType >::SizeType

◆ Third_Order_Tensor

template<class TDataType >
typedef DenseVector<Second_Order_Tensor> Kratos::GeoMechanicsMathUtilities< TDataType >::Third_Order_Tensor

◆ VectorType

template<class TDataType >
typedef Vector Kratos::GeoMechanicsMathUtilities< TDataType >::VectorType

Member Function Documentation

◆ CardanoFormula()

template<class TDataType >
static bool Kratos::GeoMechanicsMathUtilities< TDataType >::CardanoFormula ( double  a,
double  b,
double  c,
double  d,
Vector solution 
)
inlinestatic

calculates the solutions for a given cubic polynomial equation 0= a*x^3+b*x^2+c*x+d

Parameters
acoefficient
bcoefficient
ccoefficient
dcoefficient
ZeroTolnumber treated as zero
Returns
Vector of solutions WARNING only valid cubic (not quadratic, not linear, not constant) equations with three real (not complex) solutions

◆ Clipping()

template<class TDataType >
static bool Kratos::GeoMechanicsMathUtilities< TDataType >::Clipping ( std::vector< Point * > &  clipping_points,
std::vector< Point * > &  subjected_points,
std::vector< Point * > &  result_points,
double  tolerance 
)
inlinestatic

Performs clipping on the two polygons clipping_points and subjected_points (the technique used is Sutherland-Hodgman clipping) and returns the overlapping polygon result_points. The method works in 3D. Both polygons have to be convex, but they can be slightly perturbated in 3D space, this allows for performing clipping on two interpolated interfaces

Parameters
clipping_pointsvertices of clipping polygon
subjected_pointsvertices of subjected polygon
result_pointsvertices of overlapping polygon
Returns
false= no overlapping polygon, true= overlapping polygon found

◆ DeviatoricUnity() [1/2]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::DeviatoricUnity ( array_1d< double, 81 > &  Unity)
inlinestatic

Generates the fourth order deviatoric unity tensor

Parameters
Unitythe deviatoric unity (will be overwritten)

◆ DeviatoricUnity() [2/2]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::DeviatoricUnity ( std::vector< std::vector< Matrix > > &  Unity)
inlinestatic

Generates the fourth order deviatoric unity tensor

Parameters
Unitythe deviatoric unity (will be overwritten)

◆ EigenValues()

template<class TDataType >
static Vector Kratos::GeoMechanicsMathUtilities< TDataType >::EigenValues ( const Matrix A,
double  tolerance,
double  zero 
)
inlinestatic

calculates Eigenvalues of given square matrix A. The QR Algorithm with shifts is used

Parameters
Athe given square matrix the eigenvalues are to be calculated.
toleranceconvergence criteria
zeronumber treated as zero
Returns
Vector of eigenvalues WARNING only valid for 2*2 and 3*3 Matrices yet

◆ EigenValuesDirectMethod()

template<class TDataType >
static Vector Kratos::GeoMechanicsMathUtilities< TDataType >::EigenValuesDirectMethod ( const Matrix A)
inlinestatic

calculates the eigenvectiors using a direct method.

Parameters
Athe given square matrix the eigenvalues are to be calculated. WARNING only valid symmetric 3*3 Matrices

◆ EigenVectors()

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::EigenVectors ( const MatrixType A,
MatrixType vectors,
VectorType lambda,
double  zero_tolerance = 1e-9,
int  max_iterations = 10 
)
inlinestatic

calculates the eigenvectors and eigenvalues of given symmetric matrix A. The eigenvectors and eigenvalues are calculated using the iterative Gauss-Seidel-method

Parameters
Athe given symmetric matrix where the eigenvectors have to be calculated.
vectorswhere the eigenvectors will be stored
lambdawher the eigenvalues will be stored
zero_tolerancethe largest value considered to be zero
max_iterationsallowed

◆ MatrixToTensor() [1/2]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::MatrixToTensor ( MatrixType A,
array_1d< double, 81 > &  Tensor 
)
inlinestatic

Transforms a given 6*6 Matrix to a corresponing 4th order tensor

Parameters
Tensorthe given Matrix
Vectorthe Tensor

◆ MatrixToTensor() [2/2]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::MatrixToTensor ( MatrixType A,
std::vector< std::vector< Matrix > > &  Tensor 
)
inlinestatic

Transforms a given 6*6 Matrix to a corresponing 4th order tensor

Parameters
Tensorthe given Matrix
Vectorthe Tensor

◆ NormTensor()

template<class TDataType >
static double Kratos::GeoMechanicsMathUtilities< TDataType >::NormTensor ( Matrix Tensor)
inlinestatic

Builds the norm of a gibven second order tensor

Parameters
Tensorthe given second order tensor
Returns
the norm of the given tensor

◆ QRFactorization()

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::QRFactorization ( const MatrixType A,
MatrixType Q,
MatrixType R 
)
inlinestatic

calculates the QR Factorization of given square matrix A=QR. The Factorization is performed using the householder algorithm

Parameters
Athe given square matrix the factorization is to be calculated.
Qthe result matrix Q
Rthe result matrix R

◆ StrainVectorToTensor()

template<class TDataType >
static MatrixType Kratos::GeoMechanicsMathUtilities< TDataType >::StrainVectorToTensor ( const VectorType Strains)
inlinestatic

calculates Eigenvectors and the EigenValues of given square matrix A(3x3) The QR Algorithm with shifts is used

Parameters
Athe given symmetric 3x3 real matrix
Vmatrix to store the eigenvectors in rows
dmatrix to store the eigenvalues converts 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, 2*e12 ] \) fir 2D case. Hence the deviatoric components of the strain vector are divided by 2 while they are stored into the matrix
Strainsthe given strain vector
Returns
the corresponding strain tensor in matrix form

◆ TensorToMatrix() [1/3]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::TensorToMatrix ( const array_1d< double, 81 > &  Tensor,
Matrix Matrix 
)
inlinestatic

Transforms a given 4th order tensor to a corresponing 6*6 Matrix

Parameters
Tensorthe given Tensor
Vectorthe Matrix

◆ TensorToMatrix() [2/3]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::TensorToMatrix ( Fourth_Order_Tensor Tensor,
Matrix Matrix 
)
inlinestatic

◆ TensorToMatrix() [3/3]

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::TensorToMatrix ( std::vector< std::vector< Matrix > > &  Tensor,
Matrix Matrix 
)
inlinestatic

Transforms a given 4th order tensor to a corresponing 6*6 Matrix

Parameters
Tensorthe given Tensor
Vectorthe Matrix

◆ TensorToStrainVector()

template<class TDataType >
static Vector Kratos::GeoMechanicsMathUtilities< TDataType >::TensorToStrainVector ( const Matrix Tensor)
inlinestatic

◆ TensorToVector()

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::TensorToVector ( const Matrix Tensor,
Vector Vector 
)
inlinestatic

Transforms a given symmetric Tensor of second order (3*3) to a corresponing 6*1 Vector

Parameters
Tensorthe given symmetric second order tensor
Vectorthe vector

◆ VectorToTensor()

template<class TDataType >
static void Kratos::GeoMechanicsMathUtilities< TDataType >::VectorToTensor ( const Vector Stress,
Matrix Tensor 
)
inlinestatic

Transforms a given 6*1 Vector to a corresponing symmetric Tensor of second order (3*3)

Parameters
Vectorthe given vector
Tensorthe symmetric second order tensor

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