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.
Public Member Functions | List of all members
Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster > Class Template Reference

This utility calculates the exact integration necessary for the Mortar Conditions. More...

#include <exact_mortar_segmentation_utility.h>

Collaboration diagram for Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >:

Public Member Functions

bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 
Life Cycle
 ExactMortarIntegrationUtility (const SizeType IntegrationOrder=0, const double DistanceThreshold=std::numeric_limits< double >::max(), const SizeType EchoLevel=0, const double ZeroToleranceFactor=1.0, const bool ConsiderDelaunator=false)
 This is the default constructor. More...
 
virtual ~ExactMortarIntegrationUtility ()=default
 Destructor. More...
 
Operations
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, ConditionArrayListType &rConditionsPointsSlave)
 This utility computes the exact integration of the mortar condition (just the points, not the whole integration points) More...
 
bool GetExactIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, IntegrationPointsType &rIntegrationPointsSlave)
 This utility computes the exact integration of the mortar condition. More...
 
bool GetExactAreaIntegration (const GeometryType &rOriginalSlaveGeometry, const array_1d< double, 3 > &rSlaveNormal, const GeometryType &rOriginalMasterGeometry, const array_1d< double, 3 > &rMasterNormal, double &rArea)
 This utility computes the exact integration of the mortar condition and returns the area. More...
 
void GetTotalArea (const GeometryType &rOriginalSlaveGeometry, ConditionArrayListType &rConditionsPointsSlave, double &rArea)
 It returns the total area inside the integration area. More...
 
bool TestGetExactIntegration (Condition::Pointer pSlaveCond, Condition::Pointer pMasterCond, Matrix &rCustomSolution)
 This utility computes the exact integration of the mortar condition. More...
 
double TestGetExactAreaIntegration (Condition::Pointer pSlaveCond, Condition::Pointer pMasterCond)
 This utility computes the exact integration of the mortar condition and returns the area. More...
 
double TestGetExactAreaIntegration (ModelPart &rMainModelPart, Condition::Pointer pSlaveCond)
 This utility computes the exact integration of the mortar condition and returns the area. More...
 
void TestIODebug (ModelPart &rMainModelPart, const std::string IOConsidered="GiD")
 This method is used for debugging purposes. Generates a GiD mesh to check. More...
 
Access
SizeTypeGetIntegrationOrder ()
 This method gets the current mIntegrationOrder. More...
 
void SetIntegrationOrder (const SizeType IntegrationOrder)
 This method sets the current mIntegrationOrder. More...
 
doubleGetDistanceThreshold ()
 This method gets the current mDistanceThreshold. More...
 
void SetDistanceThreshold (const double DistanceThreshold)
 This method sets the current mDistanceThreshold. More...
 
SizeTypeGetEchoLevel ()
 This method gets the current mEchoLevel. More...
 
void SetEchoLevel (const SizeType EchoLevel)
 This method sets the current mEchoLevel. More...
 
doubleGetZeroToleranceFactor ()
 This method gets the current mZeroToleranceFactor. More...
 
void SetZeroToleranceFactor (const double ZeroToleranceFactor)
 This method sets the current mZeroToleranceFactor. More...
 
boolGetConsiderDelaunator ()
 This method gets the current mConsiderDelaunator. More...
 
void SetConsiderDelaunator (const bool ConsiderDelaunator)
 This method sets the current mConsiderDelaunator. More...
 

Type Definitions

typedef 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 BelongType
 The type of points belongfs to be considered. More...
 
typedef PointBelong< TNumNodes, TNumNodesMaster > PointBelongType
 The definition of the point with belonging. More...
 
typedef std::vector< array_1d< PointBelongType, TDim > > VectorArrayPointsBelong
 An array of points belong. More...
 
typedef std::vector< array_1d< PointType, TDim > > VectorArrayPoints
 A vector of points. More...
 
typedef std::conditional< TBelong, VectorArrayPointsBelong, VectorArrayPoints >::type ConditionArrayListType
 The type of array of points to be considered depending if we are interested in derivatives or not. More...
 
typedef std::vector< PointBelongTypeVectorPointsBelong
 A vector of points for derivatives. More...
 
typedef std::vector< PointTypeVectorPoints
 A vector of normal points. More...
 
typedef std::conditional< TBelong, VectorPointsBelong, VectorPoints >::type PointListType
 The type of vector of points to be considered depending if we are interested in define derivatives or not. More...
 
typedef array_1d< PointBelongType, 3 > ArrayPointsBelong
 An array of points belong. More...
 
typedef array_1d< PointType, 3 > ArrayPoints
 An array of normal points. More...
 
typedef std::conditional< TBelong, ArrayPointsBelong, ArrayPoints >::type ArrayTriangleType
 The type of arrayt of points to be used depending if we are interested in derivatives or not. More...
 
typedef Line2D2< PointLineType
 The points line geometry. More...
 
typedef Triangle3D3< PointTriangleType
 The points triangle geometry. More...
 
typedef std::conditional< TDim==2, LineType, TriangleType >::type DecompositionType
 The geometry that will be considered for decomposition. More...
 
typedef std::size_t IndexType
 The definition of the index type. More...
 
static constexpr double ZeroTolerance = std::numeric_limits<double>::epsilon()
 Definition of epsilon. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (ExactMortarIntegrationUtility)
 Pointer definition of ExactMortarIntegrationUtility. More...
 

Protected Operations

void GetIntegrationMethod ()
 Get the integration method to consider. More...
 
GeometryType::IntegrationPointsArrayType GetIntegrationTriangle ()
 Get the integration method to consider. More...
 
template<SizeType TSizeCheck = TNumNodes>
void PushBackPoints (VectorPoints &rPointList, const array_1d< bool, TSizeCheck > &rAllInside, GeometryPointType &rThisGeometry)
 This function push backs the points that are inside. More...
 
template<SizeType TSizeCheck = TNumNodes>
void PushBackPoints (VectorPointsBelong &rPointList, const array_1d< bool, TSizeCheck > &rAllInside, GeometryPointType &rThisGeometry, const PointBelongs &rThisBelongs)
 This function push backs the points that are inside. More...
 
template<SizeType TSizeCheck = TNumNodes>
void CheckInside (array_1d< bool, TSizeCheck > &rAllInside, GeometryPointType &rGeometry1, GeometryPointType &rGeometry2, const double Tolerance)
 This function checks if the points of Geometry2 are inside Geometry1. More...
 
std::vector< IndexTypeComputeAnglesIndexes (PointListType &rPointList, const array_1d< double, 3 > &rNormal) const
 This function computes the angles indexes. More...
 
void ComputeClippingIntersections (PointListType &rPointList, const GeometryPointType &rSlaveGeometry, const GeometryPointType &rMasterGeometry, const PointType &rRefCenter)
 This function computes the angles indexes. More...
 
template<class TGeometryType = GeometryType>
bool TriangleIntersections (ConditionArrayListType &rConditionsPointsSlave, PointListType &rPointList, const TGeometryType &rOriginalSlaveGeometry, const GeometryPointType &rSlaveGeometry, const GeometryPointType &rMasterGeometry, const array_1d< double, 3 > &rSlaveTangentXi, const array_1d< double, 3 > &rSlaveTangentEta, const PointType &rRefCenter, const bool IsAllInside=false)
 This function calculates the triangles intersections (this is a module, that can be used directly in the respective function) More...
 
template<SizeType TSizeCheck = TNumNodes>
static bool CheckAllInside (const array_1d< bool, TSizeCheck > &rAllInside)
 This method checks if the whole array is true. More...
 
static bool Clipping2D (PointType &rPointIntersection, const PointType &rPointOrig1, const PointType &rPointOrig2, const PointType &rPointDest1, const PointType &rPointDest2)
 This function intersects two lines in a 2D plane. More...
 
static array_1d< double, 3 > GetNormalVector2D (const array_1d< double, 3 > &rVector)
 This function calculates in 2D the normal vector to a given one. More...
 
static double AnglePoints (const PointType &rPointOrig1, const PointType &rPointOrig2, const array_1d< double, 3 > &rAxis1, const array_1d< double, 3 > &rAxis2)
 This function calculates in 2D the angle between two points. More...
 
static bool CheckPoints (const PointType &rPointOrig, const PointType &rPointDest)
 This function checks if two points are the same one. More...
 
static double FastTriagleCheck2D (const PointType &rPointOrig1, const PointType &rPointOrig2, const PointType &rPointOrig3)
 This functions calculates the determinant of a 2D triangle (using points) to check if invert the order. More...
 
static bool CheckCenterIsInside (const array_1d< double, 2 > &rAuxiliarCenterLocalCoordinates, const SizeType NumNodes=TNumNodes)
 This method checks if the center of the geometry is inside the slave geometry (to prevent convex geometries) More...
 

Detailed Description

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
class Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >

This utility calculates the exact integration necessary for the Mortar Conditions.

The utility performs a mortar segmentation in order to obtain the exact integration of the geometry intersected

Author
Vicente Mataix Ferrandiz
Template Parameters
TDimThe dimension of work
TNumNodesThe number of nodes of the slave
TBelongIf we consider belonging of nodes or not. When you do the intersections in order to get the directional derivatives you need to know where the intersections belongs to calculate the derivatives. This says between which nodes the intersection belongs
TNumNodesMasterThe number of nodes of the master

Member Typedef Documentation

◆ ArrayPoints

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef array_1d<PointType, 3> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ArrayPoints

An array of normal points.

◆ ArrayPointsBelong

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef array_1d<PointBelongType, 3> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ArrayPointsBelong

An array of points belong.

◆ ArrayTriangleType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::conditional<TBelong, ArrayPointsBelong, ArrayPoints>::type Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ArrayTriangleType

The type of arrayt of points to be used depending if we are interested in derivatives or not.

◆ BelongType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef 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 Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::BelongType

The type of points belongfs to be considered.

◆ ConditionArrayListType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::conditional<TBelong, VectorArrayPointsBelong,VectorArrayPoints>::type Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ConditionArrayListType

The type of array of points to be considered depending if we are interested in derivatives or not.

◆ DecompositionType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::conditional<TDim == 2, LineType, TriangleType>::type Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::DecompositionType

The geometry that will be considered for decomposition.

◆ IndexType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::size_t Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::IndexType

The definition of the index type.

◆ LineType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef Line2D2<Point> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::LineType

The points line geometry.

◆ PointBelongType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef PointBelong<TNumNodes, TNumNodesMaster> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::PointBelongType

The definition of the point with belonging.

◆ PointListType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::conditional<TBelong, VectorPointsBelong, VectorPoints>::type Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::PointListType

The type of vector of points to be considered depending if we are interested in define derivatives or not.

◆ TriangleType

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef Triangle3D3<Point> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::TriangleType

The points triangle geometry.

◆ VectorArrayPoints

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::vector<array_1d<PointType, TDim> > Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::VectorArrayPoints

A vector of points.

◆ VectorArrayPointsBelong

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::vector<array_1d<PointBelongType, TDim> > Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::VectorArrayPointsBelong

An array of points belong.

◆ VectorPoints

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::vector<PointType> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::VectorPoints

A vector of normal points.

◆ VectorPointsBelong

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
typedef std::vector<PointBelongType> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::VectorPointsBelong

A vector of points for derivatives.

Constructor & Destructor Documentation

◆ ExactMortarIntegrationUtility()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ExactMortarIntegrationUtility ( const SizeType  IntegrationOrder = 0,
const double  DistanceThreshold = std::numeric_limits<double>::max(),
const SizeType  EchoLevel = 0,
const double  ZeroToleranceFactor = 1.0,
const bool  ConsiderDelaunator = false 
)
inline

This is the default constructor.

Parameters
IntegrationOrderThe integration order to consider
DistanceThresholdThe maximum distance to be considered (if too far the integration will be skiped)

◆ ~ExactMortarIntegrationUtility()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
virtual Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::~ExactMortarIntegrationUtility ( )
virtualdefault

Destructor.

Member Function Documentation

◆ AnglePoints()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
static double Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::AnglePoints ( const PointType rPointOrig1,
const PointType rPointOrig2,
const array_1d< double, 3 > &  rAxis1,
const array_1d< double, 3 > &  rAxis2 
)
inlinestaticprotected

This function calculates in 2D the angle between two points.

Parameters
rPointOrig1The points from the origin geometry
rPointOrig2The points in the destination geometry
rAxis1The axis respect the angle is calculated
rAxis2The normal to the previous axis
Returns
angle The angle formed

◆ CheckAllInside()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
template<SizeType TSizeCheck = TNumNodes>
static bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::CheckAllInside ( const array_1d< bool, TSizeCheck > &  rAllInside)
inlinestaticprotected

This method checks if the whole array is true.

Parameters
rAllInsideThe nodes that are inside or not the geometry
Returns
True if all the nodes are inside, false otherwise

◆ CheckCenterIsInside()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::CheckCenterIsInside ( const array_1d< double, 2 > &  rAuxiliarCenterLocalCoordinates,
const SizeType  NumNodes = TNumNodes 
)
inlinestaticprotected

This method checks if the center of the geometry is inside the slave geometry (to prevent convex geometries)

Parameters
AuxiliarCenterLocalCoordinatesThese are the local coordinates corresponding to the center
NumNodesThe number of nodes of the geometry
Returns
True if is inside false otherwise

◆ CheckInside()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
template<SizeType TSizeCheck = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::CheckInside ( array_1d< bool, TSizeCheck > &  rAllInside,
GeometryPointType rGeometry1,
GeometryPointType rGeometry2,
const double  Tolerance 
)
inlineprotected

This function checks if the points of Geometry2 are inside Geometry1.

Parameters
rAllInsideThe nodes that are inside or not the geometry
rGeometry1The geometry where the points are checked
rGeometry2The geometry to check

◆ CheckPoints()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
static bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::CheckPoints ( const PointType rPointOrig,
const PointType rPointDest 
)
inlinestaticprotected

This function checks if two points are the same one.

Parameters
rPointOrigThe points from the origin geometry
rPointDestThe points in the destination geometry
Returns
check The check done

◆ Clipping2D()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
static bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::Clipping2D ( PointType rPointIntersection,
const PointType rPointOrig1,
const PointType rPointOrig2,
const PointType rPointDest1,
const PointType rPointDest2 
)
inlinestaticprotected

This function intersects two lines in a 2D plane.

Parameters
rPointOrig1The first point from the origin geometry
rPointOrig2The second point from the origin geometry
rPointDest1The first point in the destination geometry
rPointDest2The second point in the destination geometry
rPointIntersectionThe intersection point if there is any
Returns
True if there is a intersection point, false otherwise

◆ ComputeAnglesIndexes()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
std::vector< IndexType > Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ComputeAnglesIndexes ( PointListType rPointList,
const array_1d< double, 3 > &  rNormal 
) const
inlineprotected

This function computes the angles indexes.

Parameters
rPointListThe intersection points
rNormalThe normal vector

◆ ComputeClippingIntersections()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ComputeClippingIntersections ( PointListType rPointList,
const GeometryPointType rSlaveGeometry,
const GeometryPointType rMasterGeometry,
const PointType rRefCenter 
)
inlineprotected

This function computes the angles indexes.

Parameters
rPointListThe intersection points
rSlaveGeometryThe first (slave) geometry studied (projected)
rMasterGeometryThe second (master) geometry studied (projected)
rRefCenterThe reference point to rotate

◆ FastTriagleCheck2D()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
static double Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::FastTriagleCheck2D ( const PointType rPointOrig1,
const PointType rPointOrig2,
const PointType rPointOrig3 
)
inlinestaticprotected

This functions calculates the determinant of a 2D triangle (using points) to check if invert the order.

Parameters
rPointOrig1First point
rPointOrig2Second point
rPointOrig3Third point
Returns
The DetJ

◆ GetConsiderDelaunator()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
bool& Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetConsiderDelaunator ( )
inline

This method gets the current mConsiderDelaunator.

Returns
If considering delaunator

◆ GetDistanceThreshold()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
double& Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetDistanceThreshold ( )
inline

This method gets the current mDistanceThreshold.

Returns
The distance threshold considered

◆ GetEchoLevel()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
SizeType& Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetEchoLevel ( )
inline

This method gets the current mEchoLevel.

Returns
The echo level considered

◆ GetExactAreaIntegration()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetExactAreaIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
double rArea 
)

This utility computes the exact integration of the mortar condition and returns the area.

Parameters
rOriginalSlaveGeometryThe geometry of the slave condition
rSlaveNormalThe normal of the slave condition
rOriginalMasterGeometryThe geometry of the master condition
rMasterNormalThe normal of the master condition
rAreaThe total area integrated
Returns
True if there is a common area (the geometries intersect), false otherwise

◆ GetExactIntegration() [1/12]

bool Kratos::ExactMortarIntegrationUtility< 2, 2, false >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [2/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 3, false >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [3/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 4, false >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [4/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 3, false, 4 >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [5/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 4, false, 3 >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [6/12]

bool Kratos::ExactMortarIntegrationUtility< 2, 2, true >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [7/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 3, true >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [8/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 4, true >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [9/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 3, true, 4 >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [10/12]

bool Kratos::ExactMortarIntegrationUtility< 3, 4, true, 3 >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

◆ GetExactIntegration() [11/12]

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
ConditionArrayListType rConditionsPointsSlave 
)

This utility computes the exact integration of the mortar condition (just the points, not the whole integration points)

Parameters
rOriginalSlaveGeometryThe geometry of the slave condition
rSlaveNormalThe normal of the slave condition
rOriginalMasterGeometryThe geometry of the master condition
rMasterNormalThe normal of the master condition
rConditionsPointsSlaveThe points that perform the exact integration
Returns
True if there is a common area (the geometries intersect), false otherwise

◆ GetExactIntegration() [12/12]

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetExactIntegration ( const GeometryType rOriginalSlaveGeometry,
const array_1d< double, 3 > &  rSlaveNormal,
const GeometryType rOriginalMasterGeometry,
const array_1d< double, 3 > &  rMasterNormal,
IntegrationPointsType rIntegrationPointsSlave 
)

This utility computes the exact integration of the mortar condition.

Parameters
rOriginalSlaveGeometryThe geometry of the slave condition
rSlaveNormalThe normal of the slave condition
rOriginalMasterGeometryThe geometry of the master condition
rMasterNormalThe normal of the master condition
rIntegrationPointsSlaveThe integrations points that belong to the slave
Returns
True if there is a common area (the geometries intersect), false otherwise

◆ GetIntegrationMethod()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetIntegrationMethod
protected

Get the integration method to consider.

◆ GetIntegrationOrder()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
SizeType& Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetIntegrationOrder ( )
inline

This method gets the current mIntegrationOrder.

Returns
The integration order considered

◆ GetIntegrationTriangle()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
GeometryType::IntegrationPointsArrayType Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetIntegrationTriangle
protected

Get the integration method to consider.

◆ GetNormalVector2D()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
static array_1d<double, 3> Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetNormalVector2D ( const array_1d< double, 3 > &  rVector)
inlinestaticprotected

This function calculates in 2D the normal vector to a given one.

Parameters
rVectorThe vector to compute the normal
Returns
normal The normal vector

◆ GetTotalArea()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetTotalArea ( const GeometryType rOriginalSlaveGeometry,
ConditionArrayListType rConditionsPointsSlave,
double rArea 
)

It returns the total area inside the integration area.

Parameters
rOriginalSlaveGeometryThe geometry of the slave condition
rConditionsPointsSlaveThe points that perform the exact integration
rAreaThe total area integrated

◆ GetZeroToleranceFactor()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
double& Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::GetZeroToleranceFactor ( )
inline

This method gets the current mZeroToleranceFactor.

Returns
The zero tolerance factor considered

◆ KRATOS_CLASS_POINTER_DEFINITION()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::KRATOS_CLASS_POINTER_DEFINITION ( ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >  )

Pointer definition of ExactMortarIntegrationUtility.

◆ PushBackPoints() [1/2]

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
template<SizeType TSizeCheck = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::PushBackPoints ( VectorPoints rPointList,
const array_1d< bool, TSizeCheck > &  rAllInside,
GeometryPointType rThisGeometry 
)
inlineprotected

This function push backs the points that are inside.

Parameters
rPointListThe intersection points
rAllInsideThe nodes that are already known as inside the other geometry
rThisGeometryThe geometry considered

◆ PushBackPoints() [2/2]

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
template<SizeType TSizeCheck = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::PushBackPoints ( VectorPointsBelong rPointList,
const array_1d< bool, TSizeCheck > &  rAllInside,
GeometryPointType rThisGeometry,
const PointBelongs rThisBelongs 
)
inlineprotected

This function push backs the points that are inside.

Parameters
rPointListThe intersection points
rAllInsideThe nodes that are already known as inside the other geometry
rThisGeometryThe geometry considered
rThisBelongsThis determine where it belongs each intersection

◆ SetConsiderDelaunator()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::SetConsiderDelaunator ( const bool  ConsiderDelaunator)
inline

This method sets the current mConsiderDelaunator.

Parameters
ConsiderDelaunatorIf considering delaunator

◆ SetDistanceThreshold()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::SetDistanceThreshold ( const double  DistanceThreshold)
inline

This method sets the current mDistanceThreshold.

Parameters
DistanceThresholdThe distance threshold considered

◆ SetEchoLevel()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::SetEchoLevel ( const SizeType  EchoLevel)
inline

This method sets the current mEchoLevel.

Parameters
EchoLevelThe echo level considered

◆ SetIntegrationOrder()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::SetIntegrationOrder ( const SizeType  IntegrationOrder)
inline

This method sets the current mIntegrationOrder.

Parameters
IntegrationOrderThe integration order considered

◆ SetZeroToleranceFactor()

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::SetZeroToleranceFactor ( const double  ZeroToleranceFactor)
inline

This method sets the current mZeroToleranceFactor.

Parameters
ZeroToleranceFactorThe zero tolerance factor considered

◆ TestGetExactAreaIntegration() [1/2]

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
double Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::TestGetExactAreaIntegration ( Condition::Pointer  pSlaveCond,
Condition::Pointer  pMasterCond 
)

This utility computes the exact integration of the mortar condition and returns the area.

Parameters
pSlaveCondThe slave condition
pMasterCondThe master condition
Returns
The total area integrated

◆ TestGetExactAreaIntegration() [2/2]

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
double Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::TestGetExactAreaIntegration ( ModelPart rMainModelPart,
Condition::Pointer  pSlaveCond 
)

This utility computes the exact integration of the mortar condition and returns the area.

Parameters
rMainModelPartThe main model part
pSlaveCondThe slave condition
Returns
The total area integrated

◆ TestGetExactIntegration()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::TestGetExactIntegration ( Condition::Pointer  pSlaveCond,
Condition::Pointer  pMasterCond,
Matrix rCustomSolution 
)

This utility computes the exact integration of the mortar condition.

Parameters
pSlaveCondThe slave condition
pMasterCondThe master condition
rCustomSolutionThe matrix containing the integrations points that belong to the slave
Returns
True if there is a common area (the geometries intersect), false otherwise

◆ TestIODebug()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
void Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::TestIODebug ( ModelPart rMainModelPart,
const std::string  IOConsidered = "GiD" 
)

This method is used for debugging purposes. Generates a GiD mesh to check.

Parameters
rMainModelPartThe main model part
IOConsideredThe IO considered

◆ TriangleIntersections()

template<SizeType TDim, SizeType TNumNodes, bool TBelong, SizeType TNumNodesMaster>
template<class TGeometryType >
bool Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::TriangleIntersections ( ConditionArrayListType rConditionsPointsSlave,
PointListType rPointList,
const TGeometryType &  rOriginalSlaveGeometry,
const GeometryPointType rSlaveGeometry,
const GeometryPointType rMasterGeometry,
const array_1d< double, 3 > &  rSlaveTangentXi,
const array_1d< double, 3 > &  rSlaveTangentEta,
const PointType rRefCenter,
const bool  IsAllInside = false 
)
inlineprotected

This function calculates the triangles intersections (this is a module, that can be used directly in the respective function)

Parameters
rConditionsPointsSlaveThe final solution vector, containing all the nodes
rPointListThe intersection points
rSlaveGeometryThe first (slave) geometry studied (projected)
rMasterGeometryThe second (master) geometry studied (projected)
rSlaveTangentXiThe first vector used as base to rotate
rSlaveTangentEtaThe second vector used as base to rotate
rRefCenterThe reference point to rotate
IsAllInsideTo simplify and consider the point_list directly
Returns
If there is intersection or not (true/false)

Member Data Documentation

◆ ZeroTolerance

template<SizeType TDim, SizeType TNumNodes, bool TBelong = false, SizeType TNumNodesMaster = TNumNodes>
constexpr double Kratos::ExactMortarIntegrationUtility< TDim, TNumNodes, TBelong, TNumNodesMaster >::ZeroTolerance = std::numeric_limits<double>::epsilon()
staticconstexpr

Definition of epsilon.


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