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::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster > Class Template Reference

This utilities are used in order to compute the directional derivatives during mortar contact. More...

#include <derivatives_utilities.h>

Collaboration diagram for Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >:

Type Definitions

using IndexType = std::size_t
 The index type. More...
 
using GeometryType = Geometry< Node >
 The geometry of nodes. More...
 
using NodesArrayType = Geometry< Node >::PointsArrayType
 The array of nodes contained in a geometry. More...
 
using BelongType = typename std::conditional< TNumNodes==2, PointBelongsLine2D2N, typename std::conditional< TNumNodes==3, typename std::conditional< TNumNodesMaster==3, PointBelongsTriangle3D3N, PointBelongsTriangle3D3NQuadrilateral3D4N >::type, typename std::conditional< TNumNodesMaster==3, PointBelongsQuadrilateral3D4NTriangle3D3N, PointBelongsQuadrilateral3D4N >::type >::type >::type
 The belong type (for derivatives definition) More...
 
using PointBelongType = PointBelong< TNumNodes, TNumNodesMaster >
 The points used for derivatives definition. More...
 
using GeometryPointBelongType = Geometry< PointBelongType >
 A geometry defined by the point belongs (the points used for derivatives definition) More...
 
using ConditionArrayType = array_1d< PointBelongType, TDim >
 An array of belong point to define the geometries of belong points. More...
 
using ConditionArrayListType = typename std::vector< ConditionArrayType >
 The definition of an array pf conditions. More...
 
using LineType = Line2D2< PointType >
 The line definition. More...
 
using TriangleType = Triangle3D3< PointType >
 The triangle definition. More...
 
using DecompositionType = typename std::conditional< TDim==2, LineType, TriangleType >::type
 The geometry for decomposition (line in 2D and triangle for 3D) More...
 
using DerivativeDataType = typename std::conditional< TFrictional, DerivativeDataFrictional< TDim, TNumNodes, TNumNodesMaster >, DerivativeData< TDim, TNumNodes, TNumNodesMaster > >::type
 The derivative data type. More...
 
using GeneralVariables = MortarKinematicVariablesWithDerivatives< TDim, TNumNodes, TNumNodesMaster >
 The kinematic variables. More...
 
using AeData = DualLagrangeMultiplierOperatorsWithDerivatives< TDim, TNumNodes, TFrictional, TNumNodesMaster >
 The dual LM operators. More...
 
using MortarConditionMatrices = MortarOperatorWithDerivatives< TDim, TNumNodes, TFrictional, TNumNodesMaster >
 The mortar operators. More...
 
static constexpr double CheckThresholdCoefficient = 1.0e-12
 
static constexpr double ZeroTolerance = std::numeric_limits<double>::epsilon()
 Definition of epsilon. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (DerivativesUtilities)
 Pointer definition of DerivativesUtilities. More...
 

Operations

MatrixCalculateDeltaPosition (Matrix &rDeltaPosition, const GeometryType &rThisGeometry, const ConditionArrayType &rLocalCoordinates)
 Returns a matrix with the increment of displacements, that can be used for compute the Jacobian reference (current) configuration. More...
 
static void CalculateDeltaDetjSlave (const DecompositionType &DecompGeom, const GeneralVariables &rVariables, DerivativeDataType &rDerivativeData)
 This method is used to compute the directional derivatives of the jacobian determinant. More...
 
static array_1d< array_1d< double, 3 >, TDim *TNumNodes > GPDeltaNormalSlave (const Matrix &rJacobian, const Matrix &rDNDe)
 This method is used to compute the local increment of the normal. More...
 
static array_1d< array_1d< double, 3 >, TDim *TNumNodesMaster > GPDeltaNormalMaster (const Matrix &rJacobian, const Matrix &rDNDe)
 This method is used to compute the local increment of the normal. More...
 
static array_1d< array_1d< double, 3 >, TDim *TNumNodes > DeltaNormalCenter (const GeometryType &rThisGeometry)
 It computes the delta normal of the center of the geometry. More...
 
static void CalculateDeltaNormalSlave (array_1d< BoundedMatrix< double, TNumNodes, TDim >, TNumNodes *TDim > &rDeltaNormal, GeometryType &rThisGeometry)
 Calculates the increment of the normal and in the slave condition. More...
 
static void CalculateDeltaNormalMaster (array_1d< BoundedMatrix< double, TNumNodesMaster, TDim >, TNumNodesMaster *TDim > &rDeltaNormal, const GeometryType &rThisGeometry)
 Calculates the increment of the normal and in the master condition. More...
 
static void CalculateDeltaCellVertex (const GeneralVariables &rVariables, DerivativeDataType &rDerivativeData, const array_1d< BelongType, TDim > &rTheseBelongs, const NormalDerivativesComputation ConsiderNormalVariation, const GeometryType &rSlaveGeometry, const GeometryType &rMasterGeometry, const array_1d< double, 3 > &rNormal)
 This method is used to compute the directional derivatives of the cell vertex. More...
 
static void CalculateDeltaN1 (const GeneralVariables &rVariables, DerivativeDataType &rDerivativeData, const GeometryType &rSlaveGeometry, const GeometryType &rMasterGeometry, const array_1d< double, 3 > &rSlaveNormal, const DecompositionType &rDecompGeom, const PointType &rLocalPointDecomp, const PointType &rLocalPointParent, const NormalDerivativesComputation ConsiderNormalVariation=NormalDerivativesComputation::NO_DERIVATIVES_COMPUTATION)
 Calculates the increment of the shape functions. More...
 
static void CalculateDeltaN (const GeneralVariables &rVariables, DerivativeDataType &rDerivativeData, const GeometryType &rSlaveGeometry, const GeometryType &rMasterGeometry, const array_1d< double, 3 > &rSlaveNormal, const array_1d< double, 3 > &rMasterNormal, const DecompositionType &rDecompGeom, const PointType &rLocalPointDecomp, const PointType &rLocalPointParent, const NormalDerivativesComputation ConsiderNormalVariation=NormalDerivativesComputation::NO_DERIVATIVES_COMPUTATION, const bool DualLM=false)
 Calculates the increment of the shape functions. More...
 
static MatrixCalculateDeltaPosition (Matrix &rDeltaPosition, const GeometryType &ThisGeometry)
 Returns a matrix with the increment of displacements. More...
 
static void CalculateDeltaPosition (Vector &rDeltaPosition, const GeometryType &rSlaveGeometry, const GeometryType &rMasterGeometry, const IndexType IndexNode)
 Returns a vector with the increment of displacements. More...
 
static void CalculateDeltaPosition (Vector &rDeltaPosition, const GeometryType &rSlaveGeometry, const GeometryType &rMasterGeometry, const IndexType IndexNode, const IndexType iDoF)
 Returns a vector with the increment of displacements. More...
 
static void CalculateDeltaPosition (double &rDeltaPosition, const GeometryType &rSlaveGeometry, const GeometryType &rMasterGeometry, const IndexType IndexNode, const IndexType iDoF)
 Returns a double with the increment of displacements. More...
 
static bool CalculateAeAndDeltaAe (const GeometryType &rSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rMasterGeometry, DerivativeDataType &rDerivativeData, GeneralVariables &rVariables, const NormalDerivativesComputation ConsiderNormalVariation, ConditionArrayListType &rConditionsPointsSlave, GeometryData::IntegrationMethod ThisIntegrationMethod, const double AxiSymCoeff=1.0)
 Calculate Ae and DeltaAe matrices. More...
 

Detailed Description

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
class Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >

This utilities are used in order to compute the directional derivatives during mortar contact.

The derivatives take the same argument templates than the contact conditions

Author
Vicente Mataix Ferrandiz
Gabriel Valdes Alonzo
Template Parameters
TDimThe dimension of work
TNumNodesThe number of nodes of the slave
TNormalVariationIf the normal variation is considered
TNumNodesMasterThe number of nodes of the master

Member Typedef Documentation

◆ AeData

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::AeData = DualLagrangeMultiplierOperatorsWithDerivatives<TDim, TNumNodes, TFrictional, TNumNodesMaster>

The dual LM operators.

◆ BelongType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::BelongType = typename std::conditional<TNumNodes == 2, PointBelongsLine2D2N, typename std::conditional<TNumNodes == 3, typename std::conditional<TNumNodesMaster == 3, PointBelongsTriangle3D3N, PointBelongsTriangle3D3NQuadrilateral3D4N>::type, typename std::conditional<TNumNodesMaster == 3, PointBelongsQuadrilateral3D4NTriangle3D3N, PointBelongsQuadrilateral3D4N>::type>::type>::type

The belong type (for derivatives definition)

◆ ConditionArrayListType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::ConditionArrayListType = typename std::vector<ConditionArrayType>

The definition of an array pf conditions.

◆ ConditionArrayType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::ConditionArrayType = array_1d<PointBelongType,TDim>

An array of belong point to define the geometries of belong points.

◆ DecompositionType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::DecompositionType = typename std::conditional<TDim == 2, LineType, TriangleType >::type

The geometry for decomposition (line in 2D and triangle for 3D)

◆ DerivativeDataType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::DerivativeDataType = typename std::conditional<TFrictional, DerivativeDataFrictional<TDim, TNumNodes, TNumNodesMaster>, DerivativeData<TDim, TNumNodes, TNumNodesMaster> >::type

The derivative data type.

◆ GeneralVariables

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::GeneralVariables = MortarKinematicVariablesWithDerivatives<TDim, TNumNodes, TNumNodesMaster>

The kinematic variables.

◆ GeometryPointBelongType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::GeometryPointBelongType = Geometry<PointBelongType>

A geometry defined by the point belongs (the points used for derivatives definition)

◆ GeometryType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::GeometryType = Geometry<Node>

The geometry of nodes.

◆ IndexType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::IndexType = std::size_t

The index type.

◆ LineType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::LineType = Line2D2<PointType>

The line definition.

◆ MortarConditionMatrices

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::MortarConditionMatrices = MortarOperatorWithDerivatives<TDim, TNumNodes, TFrictional, TNumNodesMaster>

The mortar operators.

◆ NodesArrayType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::NodesArrayType = Geometry<Node>::PointsArrayType

The array of nodes contained in a geometry.

◆ PointBelongType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::PointBelongType = PointBelong<TNumNodes, TNumNodesMaster>

The points used for derivatives definition.

◆ TriangleType

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
using Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::TriangleType = Triangle3D3<PointType>

The triangle definition.

Member Function Documentation

◆ CalculateAeAndDeltaAe()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
bool Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateAeAndDeltaAe ( const GeometryType rSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rMasterGeometry,
DerivativeDataType rDerivativeData,
GeneralVariables rVariables,
const NormalDerivativesComputation  ConsiderNormalVariation,
ConditionArrayListType rConditionsPointsSlave,
GeometryData::IntegrationMethod  ThisIntegrationMethod,
const double  AxiSymCoeff = 1.0 
)
static

Calculate Ae and DeltaAe matrices.

Parameters
rSlaveGeometryThe geometry of the slave side
rSlaveNormalThe normal of the slave side
rMasterGeometryThe master side geometry
rDerivativeDataThe derivatives container
rVariablesThe kinematic variables
ConsiderNormalVariationIf consider the normal derivative
rConditionsPointsSlaveThe points that configure the exact decomposition of the geometry
ThisIntegrationMethodThe integration method considered
AxiSymCoeffThe axisymmetric coefficient

◆ CalculateDeltaCellVertex()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaCellVertex ( const GeneralVariables rVariables,
DerivativeDataType rDerivativeData,
const array_1d< BelongType, TDim > &  rTheseBelongs,
const NormalDerivativesComputation  ConsiderNormalVariation,
const GeometryType rSlaveGeometry,
const GeometryType rMasterGeometry,
const array_1d< double, 3 > &  rNormal 
)
static

This method is used to compute the directional derivatives of the cell vertex.

Parameters
rVariablesThe kinematic variables
rDerivativeDataThe derivatives container
rTheseBelongsThe belongs list used in the derivatives
ConsiderNormalVariationIf consider the normal derivative
rSlaveGeometryThe slave geometry
rMasterGeometryThe master geometry
rNormalThe normal vector of the slave geometry

The procedure will be the following in order to compute the derivative of the clipping The expression of the clipping is the following: xclipp = xs1 - num/denom * diff3 Being: diff1 = xs1 - xs2; diff2 = xe2 - xs2; diff3 = xe1 - xs1; num = (diff1 x diff2) · n0 denom = (diff3 x diff2) · n0 The derivative can be defined then as: delta_xclipp = delta_xs1 - (delta_num * denom - num * delta_denom)/delta_denom**2 * diff3 - num/ denom * delta_diff3 And here: delta_num = num · delta_n0 + n0 · (delta_diff1 x diff2 + diff1 x delta_diff2) delta_denom = denom · delta_n0 + n0 · (delta_diff3 x diff2 + diff3 x delta_diff2)

◆ CalculateDeltaDetjSlave()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaDetjSlave ( const DecompositionType DecompGeom,
const GeneralVariables rVariables,
DerivativeDataType rDerivativeData 
)
static

This method is used to compute the directional derivatives of the jacobian determinant.

Parameters
DecompGeomThe triangle used to decompose the geometry
rVariablesThe kinematic variables
rDerivativeDataThe derivatives container

◆ CalculateDeltaN()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaN ( const GeneralVariables rVariables,
DerivativeDataType rDerivativeData,
const GeometryType rSlaveGeometry,
const GeometryType rMasterGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const array_1d< double, 3 > &  rMasterNormal,
const DecompositionType rDecompGeom,
const PointType rLocalPointDecomp,
const PointType rLocalPointParent,
const NormalDerivativesComputation  ConsiderNormalVariation = NormalDerivativesComputation::NO_DERIVATIVES_COMPUTATION,
const bool  DualLM = false 
)
static

Calculates the increment of the shape functions.

Parameters
rVariablesThe kinematic variables
rDerivativeDataThe derivatives container
rSlaveGeometryThe geometry of the slave side
rMasterGeometryThe geometry of the master side
rSlaveNormalThe normal of the slave side
rMasterNormalThe normal of the master side
rDecompGeomThe triangle used to decompose the geometry
rLocalPointDecompThe local coordinates in the decomposed geometry
rLocalPointParentThe local coordinates in the slave geometry
ConsiderNormalVariationIf consider the normal derivative
DualLMIf the dual Lm formulation is considered

◆ CalculateDeltaN1()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaN1 ( const GeneralVariables rVariables,
DerivativeDataType rDerivativeData,
const GeometryType rSlaveGeometry,
const GeometryType rMasterGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const DecompositionType rDecompGeom,
const PointType rLocalPointDecomp,
const PointType rLocalPointParent,
const NormalDerivativesComputation  ConsiderNormalVariation = NormalDerivativesComputation::NO_DERIVATIVES_COMPUTATION 
)
inlinestatic

Calculates the increment of the shape functions.

Parameters
rVariablesThe kinematic variables
rDerivativeDataThe derivatives container
rSlaveGeometryThe geometry of the slave side
rMasterGeometryThe geometry of the master side
rSlaveNormalThe normal of the slave side
rDecompGeomThe triangle used to decompose the geometry
rLocalPointDecompThe local coordinates in the decomposed geometry
rLocalPointParentThe local coordinates in the slave geometry
ConsiderNormalVariationIf consider the normal derivative

◆ CalculateDeltaNormalMaster()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaNormalMaster ( array_1d< BoundedMatrix< double, TNumNodesMaster, TDim >, TNumNodesMaster *TDim > &  rDeltaNormal,
const GeometryType rThisGeometry 
)
static

Calculates the increment of the normal and in the master condition.

Parameters
rDeltaNormalThe derivative of the normal
rThisGeometryThe geometry of the master side
Note
Hardcopied for performance

◆ CalculateDeltaNormalSlave()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaNormalSlave ( array_1d< BoundedMatrix< double, TNumNodes, TDim >, TNumNodes *TDim > &  rDeltaNormal,
GeometryType rThisGeometry 
)
static

Calculates the increment of the normal and in the slave condition.

Parameters
rDeltaNormalThe derivative of the normal
rThisGeometryThe geometry of the salve side

◆ CalculateDeltaPosition() [1/5]

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaPosition ( double rDeltaPosition,
const GeometryType rSlaveGeometry,
const GeometryType rMasterGeometry,
const IndexType  IndexNode,
const IndexType  iDoF 
)
inlinestatic

Returns a double with the increment of displacements.

Parameters
rDeltaPositionThe resulting double with the increment of position
rSlaveGeometryThe slave geometry
rMasterGeometryThe master geometry
IndexNodeThe node index
iDoFThe degree of freedom index

◆ CalculateDeltaPosition() [2/5]

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
Matrix & Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaPosition ( Matrix rDeltaPosition,
const GeometryType rThisGeometry,
const ConditionArrayType rLocalCoordinates 
)

Returns a matrix with the increment of displacements, that can be used for compute the Jacobian reference (current) configuration.

Parameters
rDeltaPositionThe matrix with the increment of displacements
rThisGeometryThe geometry considered
rLocalCoordinatesThe array containing the local coordinates of the exact integration segment

◆ CalculateDeltaPosition() [3/5]

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
Matrix & Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaPosition ( Matrix rDeltaPosition,
const GeometryType ThisGeometry 
)
inlinestatic

Returns a matrix with the increment of displacements.

Parameters
rDeltaPositionThe matrix with the increment of displacements
ThisGeometryThe geometry considered

◆ CalculateDeltaPosition() [4/5]

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaPosition ( Vector rDeltaPosition,
const GeometryType rSlaveGeometry,
const GeometryType rMasterGeometry,
const IndexType  IndexNode 
)
inlinestatic

Returns a vector with the increment of displacements.

Parameters
rDeltaPositionThe resulting vector with the increment of position
rSlaveGeometryThe slave geometry
rMasterGeometryThe master geometry
IndexNodeThe node index

◆ CalculateDeltaPosition() [5/5]

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
void Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CalculateDeltaPosition ( Vector rDeltaPosition,
const GeometryType rSlaveGeometry,
const GeometryType rMasterGeometry,
const IndexType  IndexNode,
const IndexType  iDoF 
)
inlinestatic

Returns a vector with the increment of displacements.

Parameters
rDeltaPositionThe resulting vector with the increment of position
rSlaveGeometryThe slave geometry
rMasterGeometryThe master geometry
IndexNodeThe node index
iDoFThe degree of freedom index

◆ DeltaNormalCenter()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
array_1d< array_1d< double, 3 >, TDim *TNumNodes > Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::DeltaNormalCenter ( const GeometryType rThisGeometry)
static

It computes the delta normal of the center of the geometry.

Parameters
rThisGeometryThe geometry where the delta normal is computed

◆ GPDeltaNormalMaster()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
array_1d< array_1d< double, 3 >, TDim *TNumNodesMaster > Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::GPDeltaNormalMaster ( const Matrix rJacobian,
const Matrix rDNDe 
)
inlinestatic

This method is used to compute the local increment of the normal.

Parameters
rJacobianThe jacobian on the GP
rDNDeThe local gradient
Returns
The matrix containing the delta normals
Note
Not the mean, look in the contact utilities
Hardcopied for performance

◆ GPDeltaNormalSlave()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster>
array_1d< array_1d< double, 3 >, TDim *TNumNodes > Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::GPDeltaNormalSlave ( const Matrix rJacobian,
const Matrix rDNDe 
)
inlinestatic

This method is used to compute the local increment of the normal.

Parameters
rJacobianThe jacobian on the GP
rDNDeThe local gradient
Returns
The matrix containing the delta normals
Note
Not the mean, look in the contact utilities

◆ KRATOS_CLASS_POINTER_DEFINITION()

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::KRATOS_CLASS_POINTER_DEFINITION ( DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >  )

Pointer definition of DerivativesUtilities.

Member Data Documentation

◆ CheckThresholdCoefficient

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
constexpr double Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::CheckThresholdCoefficient = 1.0e-12
staticconstexpr

◆ ZeroTolerance

template<const SizeType TDim, const SizeType TNumNodes, bool TFrictional, const bool TNormalVariation, const SizeType TNumNodesMaster = TNumNodes>
constexpr double Kratos::DerivativesUtilities< TDim, TNumNodes, TFrictional, TNormalVariation, TNumNodesMaster >::ZeroTolerance = std::numeric_limits<double>::epsilon()
staticconstexpr

Definition of epsilon.


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