![]() |
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.
|
Typedefs | |
Type Definitions | |
| typedef Node | NodeType |
| typedef Point | PointType |
| typedef PointType::CoordinatesArrayType | CoordinatesArrayType |
| typedef Geometry< NodeType > | GeometryType |
| Definition of geometries. More... | |
| typedef Geometry< PointType > | GeometryPointType |
| typedef std::size_t | IndexType |
| Index type definition. More... | |
| typedef std::size_t | SizeType |
| Size type definition. More... | |
| typedef std::unordered_map< IndexType, IndexType > | IntMap |
| A map for integers. More... | |
Functions | |
| array_1d< double, 3 > | GaussPointUnitNormal (const Vector &rN, const GeometryType &rGeometry) |
| This function calculates the r_normal in a specific GP with a given shape function. More... | |
| template<> | |
| void | ResetAuxiliarValue< Variable< double > > (ModelPart &rThisModelPart) |
| template<> | |
| void | ResetAuxiliarValue< Variable< array_1d< double, 3 > > > (ModelPart &rThisModelPart) |
| template<> | |
| const std::string | GetAuxiliarVariable< Variable< double > > () |
| template<> | |
| const std::string | GetAuxiliarVariable< Variable< array_1d< double, 3 > > > () |
| template<> | |
| double | GetAuxiliarValue< Variable< double > > (NodeType &rThisNode, const std::size_t iSize) |
| template<> | |
| double | GetAuxiliarValue< Variable< array_1d< double, 3 > > > (NodeType &rThisNode, const std::size_t iSize) |
Functions | |
| bool | LengthCheck (const GeometryPointType &rGeometryLine, const double Tolerance=1.0e-6) |
| This functions checks if the length of the line is to short, with the potential of provoque ill condition in the dual LM formulation. More... | |
| bool | HeronCheck (const GeometryPointType &rGeometryTriangle) |
| This functions checks if the semiperimeter is smaller than any of the sides of the triangle. More... | |
| bool | HeronCheck (const PointType &rPointOrig1, const PointType &rPointOrig2, const PointType &rPointOrig3) |
| This functions checks if the semiperimeter is smaller than any of the sides of the triangle. More... | |
| void | RotatePoint (PointType &rPointToRotate, const PointType &rPointReferenceRotation, const array_1d< double, 3 > &rSlaveTangentXi, const array_1d< double, 3 > &rSlaveTangentEta, const bool Inversed) |
| This function rotates to align the projected points to a parallel plane to XY. More... | |
| void | ComputeNodesMeanNormalModelPart (ModelPart &rModelPart, const bool ComputeConditions=true) |
| It computes the mean of the normal in the condition in all the nodes. More... | |
| void | ComputeNodesTangentModelPart (ModelPart &rModelPart, const Variable< array_1d< double, 3 >> *pSlipVariable=NULL, const double SlipCoefficient=1.0, const bool SlipAlways=false) |
| It computes the tangent in all the nodes of the model part. More... | |
| void | ComputeNodesTangentFromNormalModelPart (ModelPart &rModelPart) |
| It computes the tangent in all the nodes of the model part from its normal. More... | |
| void | ComputeTangentsFromNormal (NodeType &rNode, const array_1d< double, 3 > &rNormal, const std::size_t Dimension=3) |
| It computes the tangent on the given node using the normal provided. More... | |
| void | ComputeTangentNodeWithLMAndSlip (NodeType &rNode, const std::size_t StepLM=0, const Variable< array_1d< double, 3 >> *pSlipVariable=NULL, const double SlipCoefficient=1.0, const std::size_t Dimension=3) |
| It computes the tangent on the given node using the LM direction and Slip direction. More... | |
| void | ComputeTangentNodeWithSlip (NodeType &rNode, const std::size_t StepLM=0, const Variable< array_1d< double, 3 >> *pSlipVariable=NULL, const double SlipCoefficient=1.0, const std::size_t Dimension=3) |
| It computes the tangent on the given node using the Slip direction. More... | |
| template<typename TType > | |
| std::vector< std::size_t > | SortIndexes (const std::vector< TType > &rThisVector) |
| This function gives you the indexes needed to order a vector. More... | |
| template<class TContainerType > | |
| void | InvertNormalForFlag (TContainerType &rContainer, const Flags Flag) |
| It inverts the order of the nodes in the conditions of a model part in order to invert the normal when certain flag is active. More... | |
| template<class TContainerType > | |
| void | InvertNormal (TContainerType &rContainer) |
| It inverts the order of the nodes in the conditions of a model part in order to invert the normal. More... | |
| template<SizeType TDim, SizeType TNumNodes> | |
| BoundedMatrix< double, TNumNodes, TDim > | GetCoordinates (const GeometryType &rGeometry, const bool Current=true, const IndexType Step=0) |
| It calculates the matrix of coordinates of a geometry. More... | |
| template<SizeType TNumNodes, SizeType TDim> | |
| BoundedMatrix< double, TNumNodes, TDim > | ComputeTangentMatrix (const GeometryType &rGeometry) |
| It calculates the matrix containing the tangent vector TANGENT_XI. More... | |
| template<SizeType TNumNodes, class TVarType = Variable<double>> | |
| array_1d< double, TNumNodes > | GetVariableVector (const GeometryType &rGeometry, const TVarType &rVariable, const IndexType Step) |
| It calculates the vector of an historical variable of a geometry. More... | |
| template<SizeType TNumNodes, class TVarType = Variable<double>> | |
| BoundedMatrix< double, TNumNodes, 1 > | GetVariableVectorMatrix (const GeometryType &rGeometry, const TVarType &rVariable, const unsigned int Step) |
| It calculates the vector of an historical variable of a geometry. More... | |
| template<SizeType TNumNodes, class TVarType = Variable<double>> | |
| array_1d< double, TNumNodes > | GetVariableVector (const GeometryType &rGeometry, const TVarType &rVariable) |
| It calculates the vector of a non-historical variable of a geometry. More... | |
| template<SizeType TNumNodes, class TVarType = Variable<double>> | |
| BoundedMatrix< double, TNumNodes, 1 > | GetVariableVectorMatrix (const GeometryType &rGeometry, const TVarType &rVariable) |
| It calculates the vector of a non-historical variable of a geometry. More... | |
| template<SizeType TDim, SizeType TNumNodes> | |
| BoundedMatrix< double, TNumNodes, TDim > | GetVariableMatrix (const GeometryType &rGeometry, const Variable< array_1d< double, 3 > > &rVariable, const unsigned int Step) |
| It calculates the matrix of a variable of a geometry. More... | |
| template<SizeType TDim, SizeType TNumNodes> | |
| BoundedMatrix< double, TNumNodes, TDim > | GetVariableMatrix (const GeometryType &rGeometry, const Variable< array_1d< double, 3 > > &rVariable) |
| It calculates the matrix of a non-historical variable of a geometry. More... | |
| template<SizeType TDim, SizeType TNumNodes> | |
| BoundedMatrix< double, TNumNodes, TDim > | GetAbsMatrix (const BoundedMatrix< double, TNumNodes, TDim > &rInputMatrix) |
| It calculates the matrix containing the absolute value of another matrix. More... | |
| template<SizeType TDim, class TVarType > | |
| unsigned int | SizeToCompute () |
| This method gives the size to be computed. More... | |
| template<class TVarType , bool THistorical> | |
| void | ResetValue (ModelPart &rThisModelPart, const TVarType &rThisVariable) |
| This method resets the value. More... | |
| template<class TVarType > | |
| void | ResetAuxiliarValue (ModelPart &rThisModelPart) |
| This method resets the auxiliar value. More... | |
| template<class TVarType > | |
| const std::string | GetAuxiliarVariable () |
| This method returns the auxiliar variable. More... | |
| template<class TVarType > | |
| double | GetAuxiliarValue (NodeType &rThisNode, const std::size_t iSize) |
| This method returns the auxiliar variable. More... | |
| template<class TVarType , bool THistorical> | |
| void | MatrixValue (const GeometryType &rThisGeometry, const TVarType &rThisVariable, Matrix &rThisValue) |
| This method adds the value. More... | |
| template<class TVarType , bool THistorical> | |
| void | AddValue (GeometryType &rThisGeometry, const TVarType &rThisVariable, const Matrix &rThisValue) |
| This method adds the value. More... | |
| template<class TVarType , bool THistorical> | |
| void | AddAreaWeightedNodalValue (NodeType &rThisNode, const TVarType &rThisVariable, const double RefArea=1.0, const double Tolerance=1.0e-4) |
| This method adds the value. More... | |
| template<class TVarType , bool THistorical> | |
| void | UpdateDatabase (ModelPart &rThisModelPart, const TVarType &rThisVariable, Vector &rDx, const std::size_t Index, IntMap &rConectivityDatabase) |
| This method updates the database in the amster side. More... | |
Definition of geometries.
| typedef std::size_t Kratos::MortarUtilities::IndexType |
Index type definition.
| typedef std::unordered_map<IndexType, IndexType> Kratos::MortarUtilities::IntMap |
A map for integers.
| typedef std::size_t Kratos::MortarUtilities::SizeType |
Size type definition.
| void Kratos::MortarUtilities::AddAreaWeightedNodalValue | ( | NodeType & | rThisNode, |
| const TVarType & | rThisVariable, | ||
| const double | RefArea = 1.0, |
||
| const double | Tolerance = 1.0e-4 |
||
| ) |
This method adds the value.
| rThisNode | The node to update |
| rThisVariable | The variable to set |
| void Kratos::MortarUtilities::AddValue | ( | GeometryType & | rThisGeometry, |
| const TVarType & | rThisVariable, | ||
| const Matrix & | rThisValue | ||
| ) |
This method adds the value.
| rThisGeometry | The geometrty to update |
| rThisVariable | The variable to set |
| rThisValue | The matrix to be updated |
| void Kratos::MortarUtilities::ComputeNodesMeanNormalModelPart | ( | ModelPart & | rModelPart, |
| const bool | ComputeConditions = true |
||
| ) |
It computes the mean of the normal in the condition in all the nodes.
| rModelPart | The model part to compute |
| ComputeConditions | If computed over conditions or elements |
| void Kratos::MortarUtilities::ComputeNodesTangentFromNormalModelPart | ( | ModelPart & | rModelPart | ) |
It computes the tangent in all the nodes of the model part from its normal.
| rModelPart | The model part to compute |
| void Kratos::MortarUtilities::ComputeNodesTangentModelPart | ( | ModelPart & | rModelPart, |
| const Variable< array_1d< double, 3 >> * | pSlipVariable = NULL, |
||
| const double | SlipCoefficient = 1.0, |
||
| const bool | SlipAlways = false |
||
| ) |
It computes the tangent in all the nodes of the model part.
| rModelPart | The model part to compute |
| pSlipVariable | The pointer to the slip variable |
| SlipCoefficient | The slip contribution |
| SlipAlways | Uses the slip even in case that LM are available |
| BoundedMatrix<double, TNumNodes, TDim> Kratos::MortarUtilities::ComputeTangentMatrix | ( | const GeometryType & | rGeometry | ) |
It calculates the matrix containing the tangent vector TANGENT_XI.
| rGeometry | The geometry to calculate |
| void Kratos::MortarUtilities::ComputeTangentNodeWithLMAndSlip | ( | NodeType & | rNode, |
| const std::size_t | StepLM = 0, |
||
| const Variable< array_1d< double, 3 >> * | pSlipVariable = NULL, |
||
| const double | SlipCoefficient = 1.0, |
||
| const std::size_t | Dimension = 3 |
||
| ) |
It computes the tangent on the given node using the LM direction and Slip direction.
| rNode | The node where to compute the tangent |
| StepLM | The considered step slip |
| pSlipVariable | The pointer to the slip variable |
| SlipCoefficient | The slip contribution |
| Dimension | The current working dimension |
| void Kratos::MortarUtilities::ComputeTangentNodeWithSlip | ( | NodeType & | rNode, |
| const std::size_t | StepLM = 0, |
||
| const Variable< array_1d< double, 3 >> * | pSlipVariable = NULL, |
||
| const double | SlipCoefficient = 1.0, |
||
| const std::size_t | Dimension = 3 |
||
| ) |
It computes the tangent on the given node using the Slip direction.
| rNode | The node where to compute the tangent |
| StepLM | The considered step slip |
| pSlipVariable | The pointer to the slip variable |
| SlipCoefficient | The slip contribution |
| Dimension | The current working dimension |
| void Kratos::MortarUtilities::ComputeTangentsFromNormal | ( | NodeType & | rNode, |
| const array_1d< double, 3 > & | rNormal, | ||
| const std::size_t | Dimension = 3 |
||
| ) |
It computes the tangent on the given node using the normal provided.
| rNode | The node where to compute the tangent |
| rNormal | The normal vector |
| Dimension | The current working dimension |
| array_1d<double,3> Kratos::MortarUtilities::GaussPointUnitNormal | ( | const Vector & | rN, |
| const GeometryType & | rGeometry | ||
| ) |
This function calculates the r_normal in a specific GP with a given shape function.
| rN | The shape function considered |
| rGeometry | The geometry of condition of interest |
| BoundedMatrix<double, TNumNodes, TDim> Kratos::MortarUtilities::GetAbsMatrix | ( | const BoundedMatrix< double, TNumNodes, TDim > & | rInputMatrix | ) |
It calculates the matrix containing the absolute value of another matrix.
| rInputMatrix | The original matrix |
| double Kratos::MortarUtilities::GetAuxiliarValue | ( | NodeType & | rThisNode, |
| const std::size_t | iSize | ||
| ) |
This method returns the auxiliar variable.
| rThisNode | Reference to the node of interest |
| iSize | The Index of the component |
| double Kratos::MortarUtilities::GetAuxiliarValue< Variable< array_1d< double, 3 > > > | ( | NodeType & | rThisNode, |
| const std::size_t | iSize | ||
| ) |
| double Kratos::MortarUtilities::GetAuxiliarValue< Variable< double > > | ( | NodeType & | rThisNode, |
| const std::size_t | iSize | ||
| ) |
| const std::string Kratos::MortarUtilities::GetAuxiliarVariable | ( | ) |
This method returns the auxiliar variable.
| const std::string Kratos::MortarUtilities::GetAuxiliarVariable< Variable< array_1d< double, 3 > > > | ( | ) |
| const std::string Kratos::MortarUtilities::GetAuxiliarVariable< Variable< double > > | ( | ) |
| BoundedMatrix<double, TNumNodes, TDim> Kratos::MortarUtilities::GetCoordinates | ( | const GeometryType & | rGeometry, |
| const bool | Current = true, |
||
| const IndexType | Step = 0 |
||
| ) |
It calculates the matrix of coordinates of a geometry.
| rGeometry | The geometry to calculate |
| Current | If we calculate the Current coordinates or the initial ones |
| Step | The time step where it is computed |
| BoundedMatrix<double, TNumNodes, TDim> Kratos::MortarUtilities::GetVariableMatrix | ( | const GeometryType & | rGeometry, |
| const Variable< array_1d< double, 3 > > & | rVariable | ||
| ) |
It calculates the matrix of a non-historical variable of a geometry.
| rGeometry | The geometry to calculate |
| rVariable | The name of the variable to calculate |
| BoundedMatrix<double, TNumNodes, TDim> Kratos::MortarUtilities::GetVariableMatrix | ( | const GeometryType & | rGeometry, |
| const Variable< array_1d< double, 3 > > & | rVariable, | ||
| const unsigned int | Step | ||
| ) |
It calculates the matrix of a variable of a geometry.
| rGeometry | The geometry to calculate |
| rVariable | The name of the variable to calculate |
| Step | The step where it is computed |
| array_1d<double, TNumNodes> Kratos::MortarUtilities::GetVariableVector | ( | const GeometryType & | rGeometry, |
| const TVarType & | rVariable | ||
| ) |
It calculates the vector of a non-historical variable of a geometry.
| rGeometry | The geometry to calculate |
| rVariable | The name of the variable to calculate |
| array_1d<double, TNumNodes> Kratos::MortarUtilities::GetVariableVector | ( | const GeometryType & | rGeometry, |
| const TVarType & | rVariable, | ||
| const IndexType | Step | ||
| ) |
It calculates the vector of an historical variable of a geometry.
| rGeometry | The geometry to calculate |
| rVariable | The name of the variable to calculate |
| Step | The step where it is computed |
| BoundedMatrix<double, TNumNodes, 1> Kratos::MortarUtilities::GetVariableVectorMatrix | ( | const GeometryType & | rGeometry, |
| const TVarType & | rVariable | ||
| ) |
It calculates the vector of a non-historical variable of a geometry.
| rGeometry | The geometry to calculate |
| rVariable | The name of the variable to calculate |
| BoundedMatrix<double, TNumNodes, 1> Kratos::MortarUtilities::GetVariableVectorMatrix | ( | const GeometryType & | rGeometry, |
| const TVarType & | rVariable, | ||
| const unsigned int | Step | ||
| ) |
It calculates the vector of an historical variable of a geometry.
| rGeometry | The geometry to calculate |
| rVariable | The name of the variable to calculate |
| Step | The step where it is computed |
| bool Kratos::MortarUtilities::HeronCheck | ( | const GeometryPointType & | rGeometryTriangle | ) |
This functions checks if the semiperimeter is smaller than any of the sides of the triangle.
| rGeometryTriangle | The triangle to be checked |
| bool Kratos::MortarUtilities::HeronCheck | ( | const PointType & | rPointOrig1, |
| const PointType & | rPointOrig2, | ||
| const PointType & | rPointOrig3 | ||
| ) |
This functions checks if the semiperimeter is smaller than any of the sides of the triangle.
| rPointOrig1 | The triangle first point |
| rPointOrig2 | The triangle second point |
| rPointOrig3 | The triangle third point |
| void Kratos::MortarUtilities::InvertNormal | ( | TContainerType & | rContainer | ) |
It inverts the order of the nodes in the conditions of a model part in order to invert the normal.
| rContainer | Reference to the objective container |
| void Kratos::MortarUtilities::InvertNormalForFlag | ( | TContainerType & | rContainer, |
| const Flags | Flag | ||
| ) |
It inverts the order of the nodes in the conditions of a model part in order to invert the normal when certain flag is active.
| rContainer | Reference to the objective container |
| Flag | The flag of the entities inverted |
| bool Kratos::MortarUtilities::LengthCheck | ( | const GeometryPointType & | rGeometryLine, |
| const double | Tolerance = 1.0e-6 |
||
| ) |
This functions checks if the length of the line is to short, with the potential of provoque ill condition in the dual LM formulation.
| rGeometryLine | The line to be checked |
| Tolerance | The threshold length |
| void Kratos::MortarUtilities::MatrixValue | ( | const GeometryType & | rThisGeometry, |
| const TVarType & | rThisVariable, | ||
| Matrix & | rThisValue | ||
| ) |
This method adds the value.
| rThisGeometry | The geometrty to update |
| rThisVariable | The variable to set |
| rThisValue | The matrix to be updated |
| void Kratos::MortarUtilities::ResetAuxiliarValue | ( | ModelPart & | rThisModelPart | ) |
This method resets the auxiliar value.
| rThisModelPart | The model part to update |
| void Kratos::MortarUtilities::ResetAuxiliarValue< Variable< array_1d< double, 3 > > > | ( | ModelPart & | rThisModelPart | ) |
| void Kratos::MortarUtilities::ResetAuxiliarValue< Variable< double > > | ( | ModelPart & | rThisModelPart | ) |
| void Kratos::MortarUtilities::ResetValue | ( | ModelPart & | rThisModelPart, |
| const TVarType & | rThisVariable | ||
| ) |
This method resets the value.
| rThisModelPart | The model part to update |
| rThisVariable | The variable to set |
| void Kratos::MortarUtilities::RotatePoint | ( | PointType & | rPointToRotate, |
| const PointType & | rPointReferenceRotation, | ||
| const array_1d< double, 3 > & | rSlaveTangentXi, | ||
| const array_1d< double, 3 > & | rSlaveTangentEta, | ||
| const bool | Inversed | ||
| ) |
This function rotates to align the projected points to a parallel plane to XY.
| rPointToRotate | The points from the origin geometry and the the point rotated |
| rPointReferenceRotation | The center point used as reference to rotate |
| rSlaveTangentXi | The first tangent vector of the slave condition |
| rSlaveTangentEta | The second tangent vector of the slave condition |
| Inversed | If we rotate to the XY or we recover from XY |
| unsigned int Kratos::MortarUtilities::SizeToCompute | ( | ) |
This method gives the size to be computed.
| std::vector<std::size_t> Kratos::MortarUtilities::SortIndexes | ( | const std::vector< TType > & | rThisVector | ) |
This function gives you the indexes needed to order a vector.
| rThisVector | The vector to order |
| void Kratos::MortarUtilities::UpdateDatabase | ( | ModelPart & | rThisModelPart, |
| const TVarType & | rThisVariable, | ||
| Vector & | rDx, | ||
| const std::size_t | Index, | ||
| IntMap & | rConectivityDatabase | ||
| ) |
This method updates the database in the amster side.
| rThisModelPart | The model part |
| rThisVariable | The variable to set |
| rDx | The vector with the increment of the value |
| Index | The index used in the case of a vector variable |
| rConectivityDatabase | The database that will be used to assemble the system |