Version: 9.12.0
SMESH_MeshAlgos Namespace Reference

Data Structures

struct  Edge
 
struct  TFreeBorderPart
 
struct  CoincidentFreeBorders
 
class  Intersector
 Cut faces of a triangular mesh. More...
 
class  Triangulate
 Divide a mesh face into triangles. More...
 
struct  PolySegment
 

Typedefs

typedef std::vector< std::vector< const SMDS_MeshElement * > > TElemGroupVector
 
typedef std::vector< std::vector< const SMDS_MeshNode * > > TNodeGroupVector
 
typedef std::vector< const SMDS_MeshNode * > TFreeBorder
 
typedef std::vector< TFreeBorderTFreeBorderVec
 
typedef std::vector< TFreeBorderPartTCoincidentGroup
 
typedef std::vector< TCoincidentGroupTCoincidentGroupVec
 
typedef std::vector< std::pair< const SMDS_MeshElement *, int > > TElemIntPairVec
 
typedef std::vector< std::pair< const SMDS_MeshNode *, int > > TNodeIntPairVec
 
typedef std::vector< PolySegmentTListOfPolySegments
 

Functions

SMESH_NodeSearcherGetNodeSearcher (SMDS_Mesh &mesh)
 Return SMESH_NodeSearcher. More...
 
SMESH_NodeSearcherGetNodeSearcher (SMDS_ElemIteratorPtr elemIt)
 Return SMESH_NodeSearcher. More...
 
SMESH_ElementSearcherGetElementSearcher (SMDS_Mesh &mesh, double tolerance=-1.)
 Return SMESH_ElementSearcher. More...
 
SMESH_ElementSearcherGetElementSearcher (SMDS_Mesh &mesh, SMDS_ElemIteratorPtr elemIt, double tolerance=-1.)
 Return SMESH_ElementSearcher acting on a sub-set of elements. More...
 
bool IsOut (const SMDS_MeshElement *element, const gp_Pnt &point, double tol)
 Return true if the point is IN or ON of the element. More...
 
double GetDistance (const SMDS_MeshElement *elem, const gp_Pnt &point, gp_XYZ *closestPnt=0)
 Return minimal distance from a point to an element. More...
 
double GetDistance (const SMDS_MeshEdge *edge, const gp_Pnt &point, gp_XYZ *closestPnt=0)
 Return minimal distance from a point to an edge. More...
 
double GetDistance (const SMDS_MeshFace *face, const gp_Pnt &point, gp_XYZ *closestPnt=0)
 Return minimal distance from a point to a face. More...
 
double GetDistance (const SMDS_MeshVolume *volume, const gp_Pnt &point, gp_XYZ *closestPnt=0)
 Return minimal distance from a point to a volume. More...
 
void GetBarycentricCoords (const gp_XY &point, const gp_XY &t0, const gp_XY &t1, const gp_XY &t2, double &bc0, double &bc1)
 Returns barycentric coordinates of a point within a triangle. More...
 
bool IntersectRayVolume (const gp_Ax1 &ray, const double rayLen, const SMDS_MeshElement *vol, double &tMin, double &tMax, int &iFacetMin, int &iFacetMax)
 Intersect volume by a ray. More...
 
const SMDS_MeshElementFindFaceInSet (const SMDS_MeshNode *n1, const SMDS_MeshNode *n2, const TIDSortedElemSet &elemSet, const TIDSortedElemSet &avoidSet, int *i1=0, int *i2=0)
 Return a face having linked nodes n1 and n2 and which is. More...
 
bool FaceNormal (const SMDS_MeshElement *F, gp_XYZ &normal, bool normalized=true)
 Calculate normal of a mesh face. More...
 
int NbCommonNodes (const SMDS_MeshElement *e1, const SMDS_MeshElement *e2)
 Return number of nodes common to two elements. More...
 
std::vector< const SMDS_MeshNode * > GetCommonNodes (const SMDS_MeshElement *e1, const SMDS_MeshElement *e2)
 Return nodes common to two elements. More...
 
bool IsOn2DBoundary (const SMDS_MeshNode *node, std::vector< const SMDS_MeshNode * > *neibors=nullptr)
 Return true if a node is on a boundary of 2D mesh. More...
 
bool IsRightOrder (const SMDS_MeshElement *face, const SMDS_MeshNode *node0, const SMDS_MeshNode *node1)
 Return true if node1 encounters first in the face and node2, after. More...
 
void Get1DBranches (SMDS_ElemIteratorPtr edgeIt, TElemGroupVector &edgeGroups, TNodeGroupVector &nodeGroups, const SMDS_MeshNode *startNode=0)
 Partition given 1D elements into groups of contiguous edges. More...
 
template<class ElemIter >
void MarkElems (ElemIter it, const bool isMarked)
 Mark elements given by SMDS_Iterator. More...
 
template<class ElemIter >
void MarkElems (ElemIter it, ElemIter end, const bool isMarked)
 Mark elements given by std iterators. More...
 
template<class ElemIter >
void MarkElemNodes (ElemIter it, const bool isMarked, const bool markElem=false)
 Mark nodes of elements given by SMDS_Iterator. More...
 
template<class ElemIter >
void MarkElemNodes (ElemIter it, ElemIter end, const bool isMarked, const bool markElem=false)
 Mark elements given by std iterators. More...
 
std::vector< EdgeFindSharpEdges (SMDS_Mesh *mesh, double angle, bool addExisting)
 Return sharp edges of faces and non-manifold ones. More...
 
std::vector< std::vector< const SMDS_MeshElement * > > SeparateFacesByEdges (SMDS_Mesh *mesh, const std::vector< Edge > &edges)
 Distribute all faces of the mesh between groups using given edges. More...
 
void FindCoincidentFreeBorders (SMDS_Mesh &mesh, double tolerance, CoincidentFreeBorders &foundFreeBordes)
 Returns TFreeBorder's coincident within the given tolerance. More...
 
void FindFreeBorders (SMDS_Mesh &mesh, TFreeBorderVec &foundFreeBordes, const bool closedOnly, bool *isManifold=0, bool *isGoodOri=0)
 Returns all or only closed TFreeBorder's. More...
 
void FillHole (const TFreeBorder &freeBorder, SMDS_Mesh &mesh, std::vector< const SMDS_MeshElement * > &newFaces)
 Fill a hole defined by a TFreeBorder with 2D elements. More...
 
void DeMerge (const SMDS_MeshElement *elem, std::vector< const SMDS_MeshNode * > &newNodes, std::vector< const SMDS_MeshNode * > &noMergeNodes)
 Find nodes whose merge makes the element invalid. More...
 
SMDS_MeshMakeOffset (SMDS_ElemIteratorPtr faceIt, SMDS_Mesh &mesh, const double offset, const bool theFixIntersections, TElemIntPairVec &new2OldFaces, TNodeIntPairVec &new2OldNodes)
 Create an offset mesh of given faces. More...
 
void MakePolyLine (SMDS_Mesh *mesh, TListOfPolySegments &segments, std::vector< const SMDS_MeshElement * > &newEdges, std::vector< const SMDS_MeshNode * > &newNodes, SMDS_MeshGroup *group=0, SMESH_ElementSearcher *searcher=0)
 Create a polyline consisting of 1D mesh elements each lying on a 2D element of the initial mesh. More...
 
std::vector< EdgeMakeSlot (SMDS_ElemIteratorPtr segmentIt, double width, SMDS_Mesh *mesh, std::vector< SMDS_MeshGroup * > &groupsToUpdate)
 Create a slot of given width around given 1D elements lying on a triangle mesh. More...
 
int MaxIndex (const gp_XYZ &x)
 Return coordinate index with maximal abs value. More...
 

Typedef Documentation

◆ TCoincidentGroup

◆ TCoincidentGroupVec

◆ TElemGroupVector

typedef std::vector< std::vector< const SMDS_MeshElement* > > SMESH_MeshAlgos::TElemGroupVector

◆ TElemIntPairVec

typedef std::vector< std::pair< const SMDS_MeshElement*, int > > SMESH_MeshAlgos::TElemIntPairVec

◆ TFreeBorder

typedef std::vector<const SMDS_MeshNode*> SMESH_MeshAlgos::TFreeBorder

◆ TFreeBorderVec

◆ TListOfPolySegments

◆ TNodeGroupVector

typedef std::vector< std::vector< const SMDS_MeshNode* > > SMESH_MeshAlgos::TNodeGroupVector

◆ TNodeIntPairVec

typedef std::vector< std::pair< const SMDS_MeshNode*, int > > SMESH_MeshAlgos::TNodeIntPairVec

Function Documentation

◆ DeMerge()

void SMESH_MeshAlgos::DeMerge ( const SMDS_MeshElement elem,
std::vector< const SMDS_MeshNode * > &  newNodes,
std::vector< const SMDS_MeshNode * > &  noMergeNodes 
)

Find nodes whose merge makes the element invalid.

(Degenerated elem is OK)

Parameters
[in]elem- the element
[in]newNodes- nodes of the element after the merge
[out]noMergeNodes- nodes to undo merge

References SMDS_MeshCell::applyInterlace(), SMDS_MeshElement::GetEntityType(), SMDS_MeshElement::GetType(), SMDS_MeshCell::interlacedSmdsOrder(), SMDS_MeshElement::IsQuadratic(), SMDS_VolumeTool::Set(), SMDSAbs_Edge, SMDSAbs_Face, and SMDSAbs_Volume.

Referenced by SMESH_MeshEditor::applyMerge().

◆ FaceNormal()

◆ FillHole()

void SMESH_MeshAlgos::FillHole ( const TFreeBorder freeBorder,
SMDS_Mesh mesh,
std::vector< const SMDS_MeshElement * > &  newFaces 
)

Fill a hole defined by a TFreeBorder with 2D elements.

Fill with 2D elements a hole defined by a TFreeBorder.

References ObjectPool< X >::getNew(), and SMDSAbs_Face.

Referenced by SMESH_MeshEditor_i::FillHole().

◆ FindCoincidentFreeBorders()

◆ FindFaceInSet()

◆ FindFreeBorders()

void SMESH_MeshAlgos::FindFreeBorders ( SMDS_Mesh mesh,
TFreeBorderVec foundFreeBordes,
const bool  closedOnly,
bool *  isManifold = 0,
bool *  isGoodOri = 0 
)

◆ FindSharpEdges()

◆ Get1DBranches()

void SMESH_MeshAlgos::Get1DBranches ( SMDS_ElemIteratorPtr  theEdgeIt,
TElemGroupVector theEdgeGroups,
TNodeGroupVector theNodeGroups,
const SMDS_MeshNode theStartNode = 0 
)

Partition given 1D elements into groups of contiguous edges.

A node where number of meeting edges != 2 is a group end. An optional startNode is used to orient groups it belongs to.

Returns
a list of edge groups and a list of corresponding node groups. If a group is closed, the first and last nodes of the group are same.

References begin(), end(), SMDS_MeshNode::GetNode(), SMDS_MeshElement::GetNode(), SMDS_MeshElement::GetType(), reverse(), and SMDSAbs_Edge.

Referenced by SMESH_MeshEditor::ExtrusionAlongTrack(), and SMESH_MeshEditor_i::Get1DBranches().

◆ GetBarycentricCoords()

void SMESH_MeshAlgos::GetBarycentricCoords ( const gp_XY &  p,
const gp_XY &  t0,
const gp_XY &  t1,
const gp_XY &  t2,
double &  bc0,
double &  bc1 
)

Returns barycentric coordinates of a point within a triangle.

A not returned bc2 = 1. - bc0 - bc1. The point lies within the triangle if ( bc0 >= 0 && bc1 >= 0 && bc0+bc1 <= 1 )

Referenced by FaceQuadStruct::findCell(), SMESH_Delaunay::FindTriangle(), FaceQuadStruct::isNear(), and SMESH_MeshAlgos::Intersector::Algo::isPointInTriangle().

◆ GetCommonNodes()

std::vector< const SMDS_MeshNode * > SMESH_MeshAlgos::GetCommonNodes ( const SMDS_MeshElement e1,
const SMDS_MeshElement e2 
)

Return nodes common to two elements.

References SMDS_MeshElement::GetNode(), SMDS_MeshElement::GetNodeIndex(), and SMDS_MeshElement::NbNodes().

◆ GetDistance() [1/4]

double SMESH_MeshAlgos::GetDistance ( const SMDS_MeshEdge edge,
const gp_Pnt &  point,
gp_XYZ *  closestPnt = 0 
)

Return minimal distance from a point to an edge.

References MESHCUT::d, SMDS_MeshCell::interlacedNodesIterator(), and SMDS_MeshCell::NbNodes().

◆ GetDistance() [2/4]

double SMESH_MeshAlgos::GetDistance ( const SMDS_MeshElement elem,
const gp_Pnt &  point,
gp_XYZ *  closestPnt = 0 
)

◆ GetDistance() [3/4]

double SMESH_MeshAlgos::GetDistance ( const SMDS_MeshFace face,
const gp_Pnt &  point,
gp_XYZ *  closestPnt = 0 
)

Return minimal distance from a point to a face.

Currently we ignore non-planarity and 2nd order of face

References GetDistance(), SMESH_MeshAlgos::Triangulate::GetTriangles(), SMDS_MeshCell::NbCornerNodes(), and SMDS_MeshCell::nodesIterator().

◆ GetDistance() [4/4]

double SMESH_MeshAlgos::GetDistance ( const SMDS_MeshVolume volume,
const gp_Pnt &  point,
gp_XYZ *  closestPnt = 0 
)

◆ GetElementSearcher() [1/2]

◆ GetElementSearcher() [2/2]

SMESH_ElementSearcher * SMESH_MeshAlgos::GetElementSearcher ( SMDS_Mesh mesh,
SMDS_ElemIteratorPtr  elemIt,
double  tolerance = -1. 
)

Return SMESH_ElementSearcher acting on a sub-set of elements.

◆ GetNodeSearcher() [1/2]

SMESH_NodeSearcher * SMESH_MeshAlgos::GetNodeSearcher ( SMDS_ElemIteratorPtr  elemIt)

◆ GetNodeSearcher() [2/2]

SMESH_NodeSearcher * SMESH_MeshAlgos::GetNodeSearcher ( SMDS_Mesh mesh)

◆ IntersectRayVolume()

bool SMESH_MeshAlgos::IntersectRayVolume ( const gp_Ax1 &  ray,
const double  rayLen,
const SMDS_MeshElement vol,
double &  tMin,
double &  tMax,
int &  iFacetMin,
int &  iFacetMax 
)

Intersect volume by a ray.

Intersect a ray with a convex volume.

Parameters
[in]ray- the ray
[in]rayLen- ray length
[in]vol- the volume
[out]tMin- return a ray parameter where the ray enters the volume
[out]tMax- return a ray parameter where the ray exit the volume
[out]iFacetMin- facet index where the ray enters the volume
[out]iFacetMax- facet index where the ray exit the volume
Returns
bool - true if the ray intersects the volume

References SMDS_VolumeTool::GetFaceNodes(), SMDS_VolumeTool::GetFaceNormal(), SMDS_VolumeTool::NbFaces(), and SMDS_VolumeTool::Set().

◆ IsOn2DBoundary()

bool SMESH_MeshAlgos::IsOn2DBoundary ( const SMDS_MeshNode theNode,
std::vector< const SMDS_MeshNode * > *  theNeibors = nullptr 
)

Return true if a node is on a boundary of 2D mesh.

Optionally returns two neighboring boundary nodes (or more in non-manifold mesh)

References SMDS_MeshNode::GetInverseElementIterator(), SMDS_MeshElement::GetNode(), SMDS_MeshElement::GetNodeIndex(), SMDS_MeshElement::NbCornerNodes(), and SMDSAbs_Face.

Referenced by SMESH_MeshEditor::RemoveNodeWithReconnection().

◆ IsOut()

◆ IsRightOrder()

bool SMESH_MeshAlgos::IsRightOrder ( const SMDS_MeshElement face,
const SMDS_MeshNode node0,
const SMDS_MeshNode node1 
)

Return true if node1 encounters first in the face and node2, after.

The nodes are supposed to be neighbor nodes in the face.

References SMDS_MeshElement::GetNodeIndex(), SMDS_MeshElement::IsMediumNode(), SMDS_MeshElement::IsQuadratic(), and SMDS_MeshElement::NbNodes().

Referenced by FindFreeBorders().

◆ MakeOffset()

SMDS_Mesh * SMESH_MeshAlgos::MakeOffset ( SMDS_ElemIteratorPtr  theFaceIt,
SMDS_Mesh theSrcMesh,
const double  theOffset,
const bool  theFixIntersections,
TElemIntPairVec theNew2OldFaces,
TNodeIntPairVec theNew2OldNodes 
)

Create an offset mesh of given faces.

Create an offsetMesh of given faces.

Parameters
[in]faceIt- the input faces
[in]theFixIntersections- to fix self intersections of the offset mesh or not
[out]new2OldFaces- history of faces
[out]new2OldNodes- history of nodes
Returns
SMDS_Mesh* - the new offset mesh, a caller should delete
Parameters
[in]faceIt- the input faces
[out]new2OldFaces- history of faces (new face -> old face ID)
[out]new2OldNodes- history of nodes (new node -> old node ID)
Returns
SMDS_Mesh* - the new offset mesh, a caller should delete

References SMDS_Mesh::AddFace(), SMDS_Mesh::AddNode(), SMDS_Mesh::AddNodeWithID(), SMDS_Mesh::AddPolygonalFace(), SMDS_Mesh::AddQuadPolygonalFace(), SMDS_MeshElement::begin_nodes(), SMDS_Mesh::ChangeElementNodes(), SMESH_MeshAlgos::Intersector::Cut(), SMDS_MeshElement::end_nodes(), FaceNormal(), SMDS_Mesh::FindNode(), GetElementSearcher(), SMDS_MeshElement::GetEntityType(), SMDS_MeshElement::GetID(), SMDS_MeshNode::GetInverseElementIterator(), SMDS_Mesh::GetMeshInfo(), SMDS_MeshElement::GetNodeIndex(), SMDS_MeshElement::GetType(), SMDS_MeshElement::isMarked(), SMESH_MeshAlgos::Intersector::MakeNewFaces(), SMDS_Mesh::MoveNode(), SMDS_MeshInfo::NbFaces(), SMDS_MeshInfo::NbTriangles(), SMESH_TNodeXYZ::Node(), SMDS_MeshElement::nodesIterator(), ORDER_QUADRATIC, SMDS_Mesh::RemoveFreeElement(), SMESH_TNodeXYZ::Set(), SMDS_MeshElement::setIsMarked(), SMDSAbs_Face, SMDSEntity_BiQuad_Quadrangle, SMDSEntity_BiQuad_Triangle, SMDSEntity_Polygon, SMDSEntity_Quad_Polygon, SMDSEntity_Quad_Quadrangle, SMDSEntity_Quad_Triangle, SMDSEntity_Quadrangle, SMDSEntity_Triangle, and SMDS_MeshNode::X().

Referenced by SMESH_MeshEditor::Offset().

◆ MakePolyLine()

void SMESH_MeshAlgos::MakePolyLine ( SMDS_Mesh mesh,
TListOfPolySegments segments,
std::vector< const SMDS_MeshElement * > &  newEdges,
std::vector< const SMDS_MeshNode * > &  newNodes,
SMDS_MeshGroup group = 0,
SMESH_ElementSearcher searcher = 0 
)

Create a polyline consisting of 1D mesh elements each lying on a 2D element of the initial mesh.

Positions of new nodes are found by cutting the mesh by the plane passing through pairs of points specified by each PolySegment structure. If there are several paths connecting a pair of points, the shortest path is selected by the module. Position of the cutting plane is defined by the two points and an optional vector lying on the plane specified by a PolySegment. By default the vector is defined by Mesh module as following. A middle point of the two given points is computed. The middle point is projected to the mesh. The vector goes from the middle point to the projection point. In case of planar mesh, the vector is normal to the mesh.

Parameters
[in,out]segments- PolySegment's defining positions of cutting planes. Return the used vector and position of the middle point.
[in]group- an optional group where created mesh segments will be added.

References SMESHUtils::Deleter< TOBJ >::_obj, SMDS_MeshGroup::Add(), SMDS_Mesh::AddEdge(), SMDS_Mesh::AddNode(), FaceNormal(), GetElementSearcher(), SMESH_ElementSearcher::GetElementsInSphere(), SMESH_MeshAlgos::PolySegment::myFace, SMESH_MeshAlgos::PolySegment::myMidProjPoint, SMESH_MeshAlgos::PolySegment::myNode1, SMESH_MeshAlgos::PolySegment::myNode2, SMESH_MeshAlgos::PolySegment::myVector, SMESH_MeshAlgos::PolySegment::myXYZ, SMDS_Mesh::nodesIterator(), SMESH_ElementSearcher::Project(), and SMDSAbs_Face.

Referenced by SMESH_MeshEditor_i::MakePolyLine(), and MakeSlot().

◆ MakeSlot()

◆ MarkElemNodes() [1/2]

template<class ElemIter >
void SMESH_MeshAlgos::MarkElemNodes ( ElemIter  it,
const bool  isMarked,
const bool  markElem = false 
)

◆ MarkElemNodes() [2/2]

template<class ElemIter >
void SMESH_MeshAlgos::MarkElemNodes ( ElemIter  it,
ElemIter  end,
const bool  isMarked,
const bool  markElem = false 
)

Mark elements given by std iterators.

References end(), and MarkElems().

◆ MarkElems() [1/2]

template<class ElemIter >
void SMESH_MeshAlgos::MarkElems ( ElemIter  it,
const bool  isMarked 
)

◆ MarkElems() [2/2]

template<class ElemIter >
void SMESH_MeshAlgos::MarkElems ( ElemIter  it,
ElemIter  end,
const bool  isMarked 
)

Mark elements given by std iterators.

References end().

◆ MaxIndex()

int SMESH_MeshAlgos::MaxIndex ( const gp_XYZ &  x)

Return coordinate index with maximal abs value.

Referenced by SMESH_MeshAlgos::Intersector::Algo::Cut(), and SMESH_MeshAlgos::Intersector::Algo::setPlaneIndices().

◆ NbCommonNodes()

int SMESH_MeshAlgos::NbCommonNodes ( const SMDS_MeshElement e1,
const SMDS_MeshElement e2 
)

◆ SeparateFacesByEdges()

std::vector< std::vector< const SMDS_MeshElement * > > SMESH_MeshAlgos::SeparateFacesByEdges ( SMDS_Mesh mesh,
const std::vector< Edge > &  edges 
)

Distribute all faces of the mesh between groups using given edges.

Distribute all faces of the mesh between groups using given edges as group boundaries.

References SMDS_MeshElement::begin_nodes(), SMDS_MeshElement::end_nodes(), SMDS_Mesh::facesIterator(), SMDS_MeshElement::isMarked(), SMDS_MeshElement::NbCornerNodes(), SMDS_Mesh::NbFaces(), and SMDS_MeshElement::setIsMarked().

Referenced by SMESH_Mesh_i::FaceGroupsSeparatedByEdges().