1 // Copyright (C) 2007-2020 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 // File : SMESH_MeshAlgos.hxx
23 // Created : Tue Apr 30 18:00:36 2013
24 // Author : Edward AGAPOV (eap)
26 // Initially this file held some low level algorithms extracted from SMESH_MeshEditor
27 // to make them accessible from Controls package
29 #include "SMESH_MeshAlgos.hxx"
31 #include "ObjectPool.hxx"
32 #include "SMDS_FaceOfNodes.hxx"
33 #include "SMDS_LinearEdge.hxx"
34 #include "SMDS_Mesh.hxx"
35 #include "SMDS_PolygonalFaceOfNodes.hxx"
36 #include "SMDS_VolumeTool.hxx"
37 #include "SMESH_OctreeNode.hxx"
39 #include <Utils_SALOME_Exception.hxx>
41 #include <GC_MakeSegment.hxx>
42 #include <GeomAPI_ExtremaCurveCurve.hxx>
43 #include <Geom_Line.hxx>
44 #include <IntAna_IntConicQuad.hxx>
45 #include <IntAna_Quadric.hxx>
48 #include <NCollection_DataMap.hxx>
53 #include <boost/container/flat_set.hpp>
55 //=======================================================================
57 * \brief Implementation of search for the node closest to point
59 //=======================================================================
61 struct SMESH_NodeSearcherImpl: public SMESH_NodeSearcher
63 //---------------------------------------------------------------------
67 SMESH_NodeSearcherImpl( const SMDS_Mesh* theMesh = 0,
68 SMDS_ElemIteratorPtr theElemIt = SMDS_ElemIteratorPtr() )
70 myMesh = ( SMDS_Mesh* ) theMesh;
72 TIDSortedNodeSet nodes;
74 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator();
76 nodes.insert( nodes.end(), nIt->next() );
80 while ( theElemIt->more() )
82 const SMDS_MeshElement* e = theElemIt->next();
83 nodes.insert( e->begin_nodes(), e->end_nodes() );
86 myOctreeNode = new SMESH_OctreeNode(nodes) ;
88 // get max size of a leaf box
89 SMESH_OctreeNode* tree = myOctreeNode;
90 while ( !tree->isLeaf() )
92 SMESH_OctreeNodeIteratorPtr cIt = tree->GetChildrenIterator();
96 myHalfLeafSize = tree->maxSize() / 2.;
99 //---------------------------------------------------------------------
101 * \brief Move node and update myOctreeNode accordingly
103 void MoveNode( const SMDS_MeshNode* node, const gp_Pnt& toPnt )
105 myOctreeNode->UpdateByMoveNode( node, toPnt );
106 myMesh->MoveNode( node, toPnt.X(), toPnt.Y(), toPnt.Z() );
109 //---------------------------------------------------------------------
113 const SMDS_MeshNode* FindClosestTo( const gp_Pnt& thePnt )
115 std::map<double, const SMDS_MeshNode*> dist2Nodes;
116 myOctreeNode->NodesAround( thePnt.Coord(), dist2Nodes, myHalfLeafSize );
117 if ( !dist2Nodes.empty() )
118 return dist2Nodes.begin()->second;
120 std::vector<const SMDS_MeshNode*> nodes;
121 //myOctreeNode->NodesAround( &tgtNode, &nodes, myHalfLeafSize );
123 double minSqDist = DBL_MAX;
124 if ( nodes.empty() ) // get all nodes of OctreeNode's closest to thePnt
126 // sort leafs by their distance from thePnt
127 typedef std::multimap< double, SMESH_OctreeNode* > TDistTreeMap;
128 TDistTreeMap treeMap;
129 std::list< SMESH_OctreeNode* > treeList;
130 std::list< SMESH_OctreeNode* >::iterator trIt;
131 treeList.push_back( myOctreeNode );
133 gp_XYZ pointNode( thePnt.X(), thePnt.Y(), thePnt.Z() );
134 bool pointInside = myOctreeNode->isInside( pointNode, myHalfLeafSize );
135 for ( trIt = treeList.begin(); trIt != treeList.end(); ++trIt)
137 SMESH_OctreeNode* tree = *trIt;
138 if ( !tree->isLeaf() ) // put children to the queue
140 if ( pointInside && !tree->isInside( pointNode, myHalfLeafSize )) continue;
141 SMESH_OctreeNodeIteratorPtr cIt = tree->GetChildrenIterator();
142 while ( cIt->more() )
143 treeList.push_back( cIt->next() );
145 else if ( tree->NbNodes() ) // put a tree to the treeMap
147 const Bnd_B3d& box = *tree->getBox();
148 double sqDist = thePnt.SquareDistance( 0.5 * ( box.CornerMin() + box.CornerMax() ));
149 treeMap.insert( std::make_pair( sqDist, tree ));
152 // find distance after which there is no sense to check tree's
153 double sqLimit = DBL_MAX;
154 TDistTreeMap::iterator sqDist_tree = treeMap.begin();
155 if ( treeMap.size() > 5 ) {
156 SMESH_OctreeNode* closestTree = sqDist_tree->second;
157 const Bnd_B3d& box = *closestTree->getBox();
158 double limit = sqrt( sqDist_tree->first ) + sqrt ( box.SquareExtent() );
159 sqLimit = limit * limit;
161 // get all nodes from trees
162 for ( ; sqDist_tree != treeMap.end(); ++sqDist_tree) {
163 if ( sqDist_tree->first > sqLimit )
165 SMESH_OctreeNode* tree = sqDist_tree->second;
166 tree->AllNodesAround( tree->GetNodeIterator()->next(), &nodes );
169 // find closest among nodes
171 const SMDS_MeshNode* closestNode = 0;
172 for ( size_t i = 0; i < nodes.size(); ++i )
174 double sqDist = thePnt.SquareDistance( SMESH_NodeXYZ( nodes[ i ]));
175 if ( minSqDist > sqDist ) {
176 closestNode = nodes[ i ];
183 //---------------------------------------------------------------------
185 * \brief Finds nodes located within a tolerance near a point
187 int FindNearPoint(const gp_Pnt& point,
188 const double tolerance,
189 std::vector< const SMDS_MeshNode* >& foundNodes)
191 myOctreeNode->NodesAround( point.Coord(), foundNodes, tolerance );
192 return foundNodes.size();
195 //---------------------------------------------------------------------
199 ~SMESH_NodeSearcherImpl() { delete myOctreeNode; }
201 //---------------------------------------------------------------------
203 * \brief Return the node tree
205 const SMESH_OctreeNode* getTree() const { return myOctreeNode; }
208 SMESH_OctreeNode* myOctreeNode;
210 double myHalfLeafSize; // max size of a leaf box
213 // ========================================================================
214 namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
216 const int MaxNbElemsInLeaf = 10; // maximal number of elements in a leaf of tree
217 const int MaxLevel = 7; // maximal tree height -> nb terminal boxes: 8^7 = 2097152
218 const double NodeRadius = 1e-9; // to enlarge bnd box of element
220 //=======================================================================
222 * \brief Octal tree of bounding boxes of elements
224 //=======================================================================
226 class ElementBndBoxTree : public SMESH_Octree
230 typedef boost::container::flat_set< const SMDS_MeshElement*, TIDCompare > TElemSeq;
232 ElementBndBoxTree(const SMDS_Mesh& mesh,
233 SMDSAbs_ElementType elemType,
234 SMDS_ElemIteratorPtr theElemIt = SMDS_ElemIteratorPtr(),
235 double tolerance = NodeRadius );
236 void getElementsNearPoint( const gp_Pnt& point, TElemSeq& foundElems );
237 void getElementsNearLine ( const gp_Ax1& line, TElemSeq& foundElems );
238 void getElementsInBox ( const Bnd_B3d& box, TElemSeq& foundElems );
239 void getElementsInSphere ( const gp_XYZ& center, const double radius, TElemSeq& foundElems );
240 ElementBndBoxTree* getLeafAtPoint( const gp_XYZ& point );
244 ElementBndBoxTree() {}
245 SMESH_Octree* newChild() const { return new ElementBndBoxTree; }
246 void buildChildrenData();
247 Bnd_B3d* buildRootBox();
249 //!< Bounding box of element
250 struct ElementBox : public Bnd_B3d
252 const SMDS_MeshElement* _element;
253 void init(const SMDS_MeshElement* elem, double tolerance);
255 std::vector< ElementBox* > _elements;
257 typedef ObjectPool< ElementBox > TElementBoxPool;
259 //!< allocator of ElementBox's and SMESH_TreeLimit
260 struct LimitAndPool : public SMESH_TreeLimit
262 TElementBoxPool _elBoPool;
263 LimitAndPool():SMESH_TreeLimit( MaxLevel, /*minSize=*/0. ) {}
265 LimitAndPool* getLimitAndPool() const
267 SMESH_TreeLimit* limitAndPool = const_cast< SMESH_TreeLimit* >( myLimit );
268 return static_cast< LimitAndPool* >( limitAndPool );
272 //================================================================================
274 * \brief ElementBndBoxTree creation
276 //================================================================================
278 ElementBndBoxTree::ElementBndBoxTree(const SMDS_Mesh& mesh,
279 SMDSAbs_ElementType elemType,
280 SMDS_ElemIteratorPtr theElemIt,
282 :SMESH_Octree( new LimitAndPool() )
284 int nbElems = mesh.GetMeshInfo().NbElements( elemType );
285 _elements.reserve( nbElems );
287 TElementBoxPool& elBoPool = getLimitAndPool()->_elBoPool;
290 if ( theElemIt && !theElemIt->more() )
291 std::cout << "WARNING: ElementBndBoxTree constructed on empty iterator!" << std::endl;
294 SMDS_ElemIteratorPtr elemIt = theElemIt ? theElemIt : mesh.elementsIterator( elemType );
295 while ( elemIt->more() )
297 ElementBox* eb = elBoPool.getNew();
298 eb->init( elemIt->next(), tolerance );
299 _elements.push_back( eb );
304 //================================================================================
306 * \brief Return the maximal box
308 //================================================================================
310 Bnd_B3d* ElementBndBoxTree::buildRootBox()
312 Bnd_B3d* box = new Bnd_B3d;
313 for ( size_t i = 0; i < _elements.size(); ++i )
314 box->Add( *_elements[i] );
318 //================================================================================
320 * \brief Redistribute element boxes among children
322 //================================================================================
324 void ElementBndBoxTree::buildChildrenData()
326 for ( size_t i = 0; i < _elements.size(); ++i )
328 for (int j = 0; j < 8; j++)
330 if ( !_elements[i]->IsOut( *myChildren[j]->getBox() ))
331 ((ElementBndBoxTree*)myChildren[j])->_elements.push_back( _elements[i]);
334 //_size = _elements.size();
335 SMESHUtils::FreeVector( _elements ); // = _elements.clear() + free memory
337 for (int j = 0; j < 8; j++)
339 ElementBndBoxTree* child = static_cast<ElementBndBoxTree*>( myChildren[j]);
340 if ((int) child->_elements.size() <= MaxNbElemsInLeaf )
341 child->myIsLeaf = true;
343 if ( child->isLeaf() && child->_elements.capacity() > child->_elements.size() )
344 SMESHUtils::CompactVector( child->_elements );
348 //================================================================================
350 * \brief Return elements which can include the point
352 //================================================================================
354 void ElementBndBoxTree::getElementsNearPoint( const gp_Pnt& point, TElemSeq& foundElems)
356 if ( getBox()->IsOut( point.XYZ() ))
361 for ( size_t i = 0; i < _elements.size(); ++i )
362 if ( !_elements[i]->IsOut( point.XYZ() ))
363 foundElems.insert( _elements[i]->_element );
367 for (int i = 0; i < 8; i++)
368 ((ElementBndBoxTree*) myChildren[i])->getElementsNearPoint( point, foundElems );
372 //================================================================================
374 * \brief Return elements which can be intersected by the line
376 //================================================================================
378 void ElementBndBoxTree::getElementsNearLine( const gp_Ax1& line, TElemSeq& foundElems )
380 if ( getBox()->IsOut( line ))
385 for ( size_t i = 0; i < _elements.size(); ++i )
386 if ( !_elements[i]->IsOut( line ) )
387 foundElems.insert( _elements[i]->_element );
391 for (int i = 0; i < 8; i++)
392 ((ElementBndBoxTree*) myChildren[i])->getElementsNearLine( line, foundElems );
396 //================================================================================
398 * \brief Return elements from leaves intersecting the sphere
400 //================================================================================
402 void ElementBndBoxTree::getElementsInSphere ( const gp_XYZ& center,
404 TElemSeq& foundElems)
406 if ( getBox()->IsOut( center, radius ))
411 for ( size_t i = 0; i < _elements.size(); ++i )
412 if ( !_elements[i]->IsOut( center, radius ))
413 foundElems.insert( _elements[i]->_element );
417 for (int i = 0; i < 8; i++)
418 ((ElementBndBoxTree*) myChildren[i])->getElementsInSphere( center, radius, foundElems );
422 //================================================================================
424 * \brief Return elements from leaves intersecting the box
426 //================================================================================
428 void ElementBndBoxTree::getElementsInBox( const Bnd_B3d& box, TElemSeq& foundElems )
430 if ( getBox()->IsOut( box ))
435 for ( size_t i = 0; i < _elements.size(); ++i )
436 if ( !_elements[i]->IsOut( box ))
437 foundElems.insert( _elements[i]->_element );
441 for (int i = 0; i < 8; i++)
442 ((ElementBndBoxTree*) myChildren[i])->getElementsInBox( box, foundElems );
446 //================================================================================
448 * \brief Return a leaf including a point
450 //================================================================================
452 ElementBndBoxTree* ElementBndBoxTree::getLeafAtPoint( const gp_XYZ& point )
454 if ( getBox()->IsOut( point ))
463 for (int i = 0; i < 8; i++)
464 if ( ElementBndBoxTree* l = ((ElementBndBoxTree*) myChildren[i])->getLeafAtPoint( point ))
470 //================================================================================
472 * \brief Return number of elements
474 //================================================================================
476 int ElementBndBoxTree::getNbElements()
481 nb = _elements.size();
485 for (int i = 0; i < 8; i++)
486 nb += ((ElementBndBoxTree*) myChildren[i])->getNbElements();
491 //================================================================================
493 * \brief Construct the element box
495 //================================================================================
497 void ElementBndBoxTree::ElementBox::init(const SMDS_MeshElement* elem, double tolerance)
500 SMDS_ElemIteratorPtr nIt = elem->nodesIterator();
501 while ( nIt->more() )
502 Add( SMESH_NodeXYZ( nIt->next() ));
503 Enlarge( tolerance );
508 //=======================================================================
510 * \brief Implementation of search for the elements by point and
511 * of classification of point in 2D mesh
513 //=======================================================================
515 SMESH_ElementSearcher::~SMESH_ElementSearcher()
519 struct SMESH_ElementSearcherImpl: public SMESH_ElementSearcher
522 SMDS_ElemIteratorPtr _meshPartIt;
523 ElementBndBoxTree* _ebbTree [SMDSAbs_NbElementTypes];
524 int _ebbTreeHeight[SMDSAbs_NbElementTypes];
525 SMESH_NodeSearcherImpl* _nodeSearcher;
526 SMDSAbs_ElementType _elementType;
528 bool _outerFacesFound;
529 std::set<const SMDS_MeshElement*> _outerFaces; // empty means "no internal faces at all"
531 SMESH_ElementSearcherImpl( SMDS_Mesh& mesh,
533 SMDS_ElemIteratorPtr elemIt=SMDS_ElemIteratorPtr())
534 : _mesh(&mesh),_meshPartIt(elemIt),_nodeSearcher(0),_tolerance(tol),_outerFacesFound(false)
536 for ( int i = 0; i < SMDSAbs_NbElementTypes; ++i )
539 _ebbTreeHeight[i] = -1;
541 _elementType = SMDSAbs_All;
543 virtual ~SMESH_ElementSearcherImpl()
545 for ( int i = 0; i < SMDSAbs_NbElementTypes; ++i )
547 delete _ebbTree[i]; _ebbTree[i] = NULL;
549 if ( _nodeSearcher ) delete _nodeSearcher;
552 virtual int FindElementsByPoint(const gp_Pnt& point,
553 SMDSAbs_ElementType type,
554 std::vector< const SMDS_MeshElement* >& foundElements);
555 virtual TopAbs_State GetPointState(const gp_Pnt& point);
556 virtual const SMDS_MeshElement* FindClosestTo( const gp_Pnt& point,
557 SMDSAbs_ElementType type );
559 virtual void GetElementsNearLine( const gp_Ax1& line,
560 SMDSAbs_ElementType type,
561 std::vector< const SMDS_MeshElement* >& foundElems);
562 virtual void GetElementsInSphere( const gp_XYZ& center,
564 SMDSAbs_ElementType type,
565 std::vector< const SMDS_MeshElement* >& foundElems);
566 virtual void GetElementsInBox( const Bnd_B3d& box,
567 SMDSAbs_ElementType type,
568 std::vector< const SMDS_MeshElement* >& foundElems);
569 virtual gp_XYZ Project(const gp_Pnt& point,
570 SMDSAbs_ElementType type,
571 const SMDS_MeshElement** closestElem);
572 double getTolerance();
573 bool getIntersParamOnLine(const gp_Lin& line, const SMDS_MeshElement* face,
574 const double tolerance, double & param);
575 void findOuterBoundary(const SMDS_MeshElement* anyOuterFace);
576 bool isOuterBoundary(const SMDS_MeshElement* face) const
578 return _outerFaces.empty() || _outerFaces.count(face);
582 if ( _ebbTreeHeight[ _elementType ] < 0 )
583 _ebbTreeHeight[ _elementType ] = _ebbTree[ _elementType ]->getHeight();
584 return _ebbTreeHeight[ _elementType ];
587 struct TInters //!< data of intersection of the line and the mesh face (used in GetPointState())
589 const SMDS_MeshElement* _face;
591 bool _coincides; //!< the line lays in face plane
592 TInters(const SMDS_MeshElement* face, const gp_Vec& faceNorm, bool coinc=false)
593 : _face(face), _faceNorm( faceNorm ), _coincides( coinc ) {}
595 struct TFaceLink //!< link and faces sharing it (used in findOuterBoundary())
598 TIDSortedElemSet _faces;
599 TFaceLink( const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshElement* face)
600 : _link( n1, n2 ), _faces( &face, &face + 1) {}
604 ostream& operator<< (ostream& out, const SMESH_ElementSearcherImpl::TInters& i)
606 return out << "TInters(face=" << ( i._face ? i._face->GetID() : 0)
607 << ", _coincides="<<i._coincides << ")";
610 //=======================================================================
612 * \brief define tolerance for search
614 //=======================================================================
616 double SMESH_ElementSearcherImpl::getTolerance()
618 if ( _tolerance < 0 )
620 const SMDS_MeshInfo& meshInfo = _mesh->GetMeshInfo();
623 if ( _nodeSearcher && meshInfo.NbNodes() > 1 )
625 double boxSize = _nodeSearcher->getTree()->maxSize();
626 _tolerance = 1e-8 * boxSize/* / meshInfo.NbNodes()*/;
628 else if ( _ebbTree[_elementType] && meshInfo.NbElements() > 0 )
630 double boxSize = _ebbTree[_elementType]->maxSize();
631 _tolerance = 1e-8 * boxSize/* / meshInfo.NbElements()*/;
633 if ( _tolerance == 0 )
635 // define tolerance by size of a most complex element
636 int complexType = SMDSAbs_Volume;
637 while ( complexType > SMDSAbs_All &&
638 meshInfo.NbElements( SMDSAbs_ElementType( complexType )) < 1 )
640 if ( complexType == SMDSAbs_All ) return 0; // empty mesh
642 if ( complexType == int( SMDSAbs_Node ))
644 SMDS_NodeIteratorPtr nodeIt = _mesh->nodesIterator();
646 if ( meshInfo.NbNodes() > 2 )
647 elemSize = SMESH_TNodeXYZ( nodeIt->next() ).Distance( nodeIt->next() );
651 SMDS_ElemIteratorPtr elemIt = _mesh->elementsIterator( SMDSAbs_ElementType( complexType ));
652 const SMDS_MeshElement* elem = elemIt->next();
653 SMDS_ElemIteratorPtr nodeIt = elem->nodesIterator();
654 SMESH_TNodeXYZ n1( nodeIt->next() );
656 while ( nodeIt->more() )
658 double dist = n1.Distance( static_cast<const SMDS_MeshNode*>( nodeIt->next() ));
659 elemSize = std::max( dist, elemSize );
662 _tolerance = 1e-4 * elemSize;
668 //================================================================================
670 * \brief Find intersection of the line and an edge of face and return parameter on line
672 //================================================================================
674 bool SMESH_ElementSearcherImpl::getIntersParamOnLine(const gp_Lin& line,
675 const SMDS_MeshElement* face,
682 GeomAPI_ExtremaCurveCurve anExtCC;
683 Handle(Geom_Curve) lineCurve = new Geom_Line( line );
685 int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
686 for ( int i = 0; i < nbNodes && nbInts < 2; ++i )
688 GC_MakeSegment edge( SMESH_TNodeXYZ( face->GetNode( i )),
689 SMESH_TNodeXYZ( face->GetNode( (i+1)%nbNodes) ));
690 anExtCC.Init( lineCurve, edge.Value() );
691 if ( anExtCC.NbExtrema() > 0 && anExtCC.LowerDistance() <= tol)
693 Standard_Real pl, pe;
694 anExtCC.LowerDistanceParameters( pl, pe );
700 if ( nbInts > 0 ) param /= nbInts;
703 //================================================================================
705 * \brief Find all faces belonging to the outer boundary of mesh
707 //================================================================================
709 void SMESH_ElementSearcherImpl::findOuterBoundary(const SMDS_MeshElement* outerFace)
711 if ( _outerFacesFound ) return;
713 // Collect all outer faces by passing from one outer face to another via their links
714 // and BTW find out if there are internal faces at all.
716 // checked links and links where outer boundary meets internal one
717 std::set< SMESH_TLink > visitedLinks, seamLinks;
719 // links to treat with already visited faces sharing them
720 std::list < TFaceLink > startLinks;
722 // load startLinks with the first outerFace
723 startLinks.push_back( TFaceLink( outerFace->GetNode(0), outerFace->GetNode(1), outerFace));
724 _outerFaces.insert( outerFace );
726 TIDSortedElemSet emptySet;
727 while ( !startLinks.empty() )
729 const SMESH_TLink& link = startLinks.front()._link;
730 TIDSortedElemSet& faces = startLinks.front()._faces;
732 outerFace = *faces.begin();
733 // find other faces sharing the link
734 const SMDS_MeshElement* f;
735 while (( f = SMESH_MeshAlgos::FindFaceInSet(link.node1(), link.node2(), emptySet, faces )))
738 // select another outer face among the found
739 const SMDS_MeshElement* outerFace2 = 0;
740 if ( faces.size() == 2 )
742 outerFace2 = (outerFace == *faces.begin() ? *faces.rbegin() : *faces.begin());
744 else if ( faces.size() > 2 )
746 seamLinks.insert( link );
748 // link direction within the outerFace
749 gp_Vec n1n2( SMESH_TNodeXYZ( link.node1()),
750 SMESH_TNodeXYZ( link.node2()));
751 int i1 = outerFace->GetNodeIndex( link.node1() );
752 int i2 = outerFace->GetNodeIndex( link.node2() );
753 bool rev = ( abs(i2-i1) == 1 ? i1 > i2 : i2 > i1 );
754 if ( rev ) n1n2.Reverse();
756 gp_XYZ ofNorm, fNorm;
757 if ( SMESH_MeshAlgos::FaceNormal( outerFace, ofNorm, /*normalized=*/false ))
759 // direction from the link inside outerFace
760 gp_Vec dirInOF = gp_Vec( ofNorm ) ^ n1n2;
761 // sort all other faces by angle with the dirInOF
762 std::map< double, const SMDS_MeshElement* > angle2Face;
763 std::set< const SMDS_MeshElement*, TIDCompare >::const_iterator face = faces.begin();
764 for ( ; face != faces.end(); ++face )
766 if ( *face == outerFace ) continue;
767 if ( !SMESH_MeshAlgos::FaceNormal( *face, fNorm, /*normalized=*/false ))
769 gp_Vec dirInF = gp_Vec( fNorm ) ^ n1n2;
770 double angle = dirInOF.AngleWithRef( dirInF, n1n2 );
771 if ( angle < 0 ) angle += 2. * M_PI;
772 angle2Face.insert( std::make_pair( angle, *face ));
774 if ( !angle2Face.empty() )
775 outerFace2 = angle2Face.begin()->second;
778 // store the found outer face and add its links to continue searching from
781 _outerFaces.insert( outerFace2 );
782 int nbNodes = outerFace2->NbCornerNodes();
783 for ( int i = 0; i < nbNodes; ++i )
785 SMESH_TLink link2( outerFace2->GetNode(i), outerFace2->GetNode((i+1)%nbNodes));
786 if ( visitedLinks.insert( link2 ).second )
787 startLinks.push_back( TFaceLink( link2.node1(), link2.node2(), outerFace2 ));
790 startLinks.pop_front();
792 _outerFacesFound = true;
794 if ( !seamLinks.empty() )
796 // There are internal boundaries touching the outher one,
797 // find all faces of internal boundaries in order to find
798 // faces of boundaries of holes, if any.
807 //=======================================================================
809 * \brief Find elements of given type where the given point is IN or ON.
810 * Returns nb of found elements and elements them-selves.
812 * 'ALL' type means elements of any type excluding nodes, balls and 0D elements
814 //=======================================================================
816 int SMESH_ElementSearcherImpl::
817 FindElementsByPoint(const gp_Pnt& point,
818 SMDSAbs_ElementType type,
819 std::vector< const SMDS_MeshElement* >& foundElements)
821 foundElements.clear();
824 double tolerance = getTolerance();
826 // =================================================================================
827 if ( type == SMDSAbs_Node || type == SMDSAbs_0DElement || type == SMDSAbs_Ball)
829 if ( !_nodeSearcher )
832 _nodeSearcher = new SMESH_NodeSearcherImpl( 0, _meshPartIt );
834 _nodeSearcher = new SMESH_NodeSearcherImpl( _mesh );
836 std::vector< const SMDS_MeshNode* > foundNodes;
837 _nodeSearcher->FindNearPoint( point, tolerance, foundNodes );
839 if ( type == SMDSAbs_Node )
841 foundElements.assign( foundNodes.begin(), foundNodes.end() );
845 for ( size_t i = 0; i < foundNodes.size(); ++i )
847 SMDS_ElemIteratorPtr elemIt = foundNodes[i]->GetInverseElementIterator( type );
848 while ( elemIt->more() )
849 foundElements.push_back( elemIt->next() );
853 // =================================================================================
854 else // elements more complex than 0D
856 if ( !_ebbTree[type] )
858 _ebbTree[_elementType] = new ElementBndBoxTree( *_mesh, type, _meshPartIt, tolerance );
860 ElementBndBoxTree::TElemSeq suspectElems;
861 _ebbTree[ type ]->getElementsNearPoint( point, suspectElems );
862 ElementBndBoxTree::TElemSeq::iterator elem = suspectElems.begin();
863 for ( ; elem != suspectElems.end(); ++elem )
864 if ( !SMESH_MeshAlgos::IsOut( *elem, point, tolerance ))
865 foundElements.push_back( *elem );
867 return foundElements.size();
870 //=======================================================================
872 * \brief Find an element of given type most close to the given point
874 * WARNING: Only edge, face and volume search is implemented so far
876 //=======================================================================
878 const SMDS_MeshElement*
879 SMESH_ElementSearcherImpl::FindClosestTo( const gp_Pnt& point,
880 SMDSAbs_ElementType type )
882 const SMDS_MeshElement* closestElem = 0;
885 if ( type == SMDSAbs_Face ||
886 type == SMDSAbs_Volume ||
887 type == SMDSAbs_Edge )
889 ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
891 ebbTree = new ElementBndBoxTree( *_mesh, type, _meshPartIt );
893 ElementBndBoxTree::TElemSeq suspectElems;
894 ebbTree->getElementsNearPoint( point, suspectElems );
896 if ( suspectElems.empty() && ebbTree->maxSize() > 0 )
898 gp_Pnt boxCenter = 0.5 * ( ebbTree->getBox()->CornerMin() +
899 ebbTree->getBox()->CornerMax() );
901 if ( ebbTree->getBox()->IsOut( point.XYZ() ))
902 radius = point.Distance( boxCenter ) - 0.5 * ebbTree->maxSize();
904 radius = ebbTree->maxSize() / pow( 2., getTreeHeight()) / 2;
905 while ( suspectElems.empty() && radius < 1e100 )
907 ebbTree->getElementsInSphere( point.XYZ(), radius, suspectElems );
911 double minDist = std::numeric_limits<double>::max();
912 std::multimap< double, const SMDS_MeshElement* > dist2face;
913 ElementBndBoxTree::TElemSeq::iterator elem = suspectElems.begin();
914 for ( ; elem != suspectElems.end(); ++elem )
916 double dist = SMESH_MeshAlgos::GetDistance( *elem, point );
917 if ( dist < minDist + 1e-10)
920 dist2face.insert( dist2face.begin(), std::make_pair( dist, *elem ));
923 if ( !dist2face.empty() )
925 std::multimap< double, const SMDS_MeshElement* >::iterator d2f = dist2face.begin();
926 closestElem = d2f->second;
927 // if there are several elements at the same distance, select one
928 // with GC closest to the point
929 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
930 double minDistToGC = 0;
931 for ( ++d2f; d2f != dist2face.end() && fabs( d2f->first - minDist ) < 1e-10; ++d2f )
933 if ( minDistToGC == 0 )
936 gc = accumulate( TXyzIterator(closestElem->nodesIterator()),
937 TXyzIterator(), gc ) / closestElem->NbNodes();
938 minDistToGC = point.SquareDistance( gc );
941 gc = accumulate( TXyzIterator( d2f->second->nodesIterator()),
942 TXyzIterator(), gc ) / d2f->second->NbNodes();
943 double d = point.SquareDistance( gc );
944 if ( d < minDistToGC )
947 closestElem = d2f->second;
950 // cout << "FindClosestTo( " <<point.X()<<", "<<point.Y()<<", "<<point.Z()<<" ) FACE "
951 // <<closestElem->GetID() << " DIST " << minDist << endl;
956 // NOT IMPLEMENTED SO FAR
962 //================================================================================
964 * \brief Classify the given point in the closed 2D mesh
966 //================================================================================
968 TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
970 _elementType = SMDSAbs_Face;
972 double tolerance = getTolerance();
974 ElementBndBoxTree*& ebbTree = _ebbTree[ SMDSAbs_Face ];
976 ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
978 // Algo: analyse transition of a line starting at the point through mesh boundary;
979 // try three lines parallel to axis of the coordinate system and perform rough
980 // analysis. If solution is not clear perform thorough analysis.
982 const int nbAxes = 3;
983 gp_Dir axisDir[ nbAxes ] = { gp::DX(), gp::DY(), gp::DZ() };
984 std::map< double, TInters > paramOnLine2TInters[ nbAxes ];
985 std::list< TInters > tangentInters[ nbAxes ]; // of faces whose plane includes the line
986 std::multimap< int, int > nbInt2Axis; // to find the simplest case
987 for ( int axis = 0; axis < nbAxes; ++axis )
989 gp_Ax1 lineAxis( point, axisDir[axis]);
990 gp_Lin line ( lineAxis );
992 ElementBndBoxTree::TElemSeq suspectFaces; // faces possibly intersecting the line
993 ebbTree->getElementsNearLine( lineAxis, suspectFaces );
995 // Intersect faces with the line
997 std::map< double, TInters > & u2inters = paramOnLine2TInters[ axis ];
998 ElementBndBoxTree::TElemSeq::iterator face = suspectFaces.begin();
999 for ( ; face != suspectFaces.end(); ++face )
1003 if ( !SMESH_MeshAlgos::FaceNormal( *face, fNorm, /*normalized=*/false)) continue;
1004 gp_Pln facePlane( SMESH_TNodeXYZ( (*face)->GetNode(0)), fNorm );
1006 // perform intersection
1007 IntAna_IntConicQuad intersection( line, IntAna_Quadric( facePlane ));
1008 if ( !intersection.IsDone() )
1010 if ( intersection.IsInQuadric() )
1012 tangentInters[ axis ].push_back( TInters( *face, fNorm, true ));
1014 else if ( ! intersection.IsParallel() && intersection.NbPoints() > 0 )
1016 double tol = 1e-4 * Sqrt( fNorm.Modulus() );
1017 gp_Pnt intersectionPoint = intersection.Point(1);
1018 if ( !SMESH_MeshAlgos::IsOut( *face, intersectionPoint, tol ))
1019 u2inters.insert( std::make_pair( intersection.ParamOnConic(1), TInters( *face, fNorm )));
1022 // Analyse intersections roughly
1024 int nbInter = u2inters.size();
1028 double f = u2inters.begin()->first, l = u2inters.rbegin()->first;
1029 if ( nbInter == 1 ) // not closed mesh
1030 return fabs( f ) < tolerance ? TopAbs_ON : TopAbs_UNKNOWN;
1032 if ( fabs( f ) < tolerance || fabs( l ) < tolerance )
1035 if ( (f<0) == (l<0) )
1038 int nbIntBeforePoint = std::distance( u2inters.begin(), u2inters.lower_bound(0));
1039 int nbIntAfterPoint = nbInter - nbIntBeforePoint;
1040 if ( nbIntBeforePoint == 1 || nbIntAfterPoint == 1 )
1043 nbInt2Axis.insert( std::make_pair( std::min( nbIntBeforePoint, nbIntAfterPoint ), axis ));
1045 if ( _outerFacesFound ) break; // pass to thorough analysis
1047 } // three attempts - loop on CS axes
1049 // Analyse intersections thoroughly.
1050 // We make two loops maximum, on the first one we only exclude touching intersections,
1051 // on the second, if situation is still unclear, we gather and use information on
1052 // position of faces (internal or outer). If faces position is already gathered,
1053 // we make the second loop right away.
1055 for ( int hasPositionInfo = _outerFacesFound; hasPositionInfo < 2; ++hasPositionInfo )
1057 std::multimap< int, int >::const_iterator nb_axis = nbInt2Axis.begin();
1058 for ( ; nb_axis != nbInt2Axis.end(); ++nb_axis )
1060 int axis = nb_axis->second;
1061 std::map< double, TInters > & u2inters = paramOnLine2TInters[ axis ];
1063 gp_Ax1 lineAxis( point, axisDir[axis]);
1064 gp_Lin line ( lineAxis );
1066 // add tangent intersections to u2inters
1068 std::list< TInters >::const_iterator tgtInt = tangentInters[ axis ].begin();
1069 for ( ; tgtInt != tangentInters[ axis ].end(); ++tgtInt )
1070 if ( getIntersParamOnLine( line, tgtInt->_face, tolerance, param ))
1071 u2inters.insert( std::make_pair( param, *tgtInt ));
1072 tangentInters[ axis ].clear();
1074 // Count intersections before and after the point excluding touching ones.
1075 // If hasPositionInfo we count intersections of outer boundary only
1077 int nbIntBeforePoint = 0, nbIntAfterPoint = 0;
1078 double f = std::numeric_limits<double>::max(), l = -std::numeric_limits<double>::max();
1079 std::map< double, TInters >::iterator u_int1 = u2inters.begin(), u_int2 = u_int1;
1080 bool ok = ! u_int1->second._coincides;
1081 while ( ok && u_int1 != u2inters.end() )
1083 double u = u_int1->first;
1084 bool touchingInt = false;
1085 if ( ++u_int2 != u2inters.end() )
1087 // skip intersections at the same point (if the line passes through edge or node)
1089 while ( u_int2 != u2inters.end() && fabs( u_int2->first - u ) < tolerance )
1095 // skip tangent intersections
1097 if ( u_int2 != u2inters.end() )
1099 const SMDS_MeshElement* prevFace = u_int1->second._face;
1100 while ( ok && u_int2->second._coincides )
1102 if ( SMESH_MeshAlgos::NbCommonNodes(prevFace , u_int2->second._face) == 0 )
1108 ok = ( u_int2 != u2inters.end() );
1114 // skip intersections at the same point after tangent intersections
1117 double u2 = u_int2->first;
1119 while ( u_int2 != u2inters.end() && fabs( u_int2->first - u2 ) < tolerance )
1125 // decide if we skipped a touching intersection
1126 if ( nbSamePnt + nbTgt > 0 )
1128 double minDot = std::numeric_limits<double>::max(), maxDot = -minDot;
1129 std::map< double, TInters >::iterator u_int = u_int1;
1130 for ( ; u_int != u_int2; ++u_int )
1132 if ( u_int->second._coincides ) continue;
1133 double dot = u_int->second._faceNorm * line.Direction();
1134 if ( dot > maxDot ) maxDot = dot;
1135 if ( dot < minDot ) minDot = dot;
1137 touchingInt = ( minDot*maxDot < 0 );
1142 if ( !hasPositionInfo || isOuterBoundary( u_int1->second._face ))
1153 u_int1 = u_int2; // to next intersection
1155 } // loop on intersections with one line
1159 if ( fabs( f ) < tolerance || fabs( l ) < tolerance )
1162 if ( nbIntBeforePoint == 0 || nbIntAfterPoint == 0)
1165 if ( nbIntBeforePoint + nbIntAfterPoint == 1 ) // not closed mesh
1166 return fabs( f ) < tolerance ? TopAbs_ON : TopAbs_UNKNOWN;
1168 if ( nbIntBeforePoint == 1 || nbIntAfterPoint == 1 )
1171 if ( (f<0) == (l<0) )
1174 if ( hasPositionInfo )
1175 return nbIntBeforePoint % 2 ? TopAbs_IN : TopAbs_OUT;
1177 } // loop on intersections of the tree lines - thorough analysis
1179 if ( !hasPositionInfo )
1181 // gather info on faces position - is face in the outer boundary or not
1182 std::map< double, TInters > & u2inters = paramOnLine2TInters[ 0 ];
1183 findOuterBoundary( u2inters.begin()->second._face );
1186 } // two attempts - with and w/o faces position info in the mesh
1188 return TopAbs_UNKNOWN;
1191 //=======================================================================
1193 * \brief Return elements possibly intersecting the line
1195 //=======================================================================
1197 void SMESH_ElementSearcherImpl::
1198 GetElementsNearLine( const gp_Ax1& line,
1199 SMDSAbs_ElementType type,
1200 std::vector< const SMDS_MeshElement* >& foundElems)
1202 _elementType = type;
1203 ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
1205 ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
1207 ElementBndBoxTree::TElemSeq elems;
1208 ebbTree->getElementsNearLine( line, elems );
1210 foundElems.insert( foundElems.end(), elems.begin(), elems.end() );
1213 //=======================================================================
1215 * Return elements whose bounding box intersects a sphere
1217 //=======================================================================
1219 void SMESH_ElementSearcherImpl::
1220 GetElementsInSphere( const gp_XYZ& center,
1221 const double radius,
1222 SMDSAbs_ElementType type,
1223 std::vector< const SMDS_MeshElement* >& foundElems)
1225 _elementType = type;
1226 ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
1228 ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
1230 ElementBndBoxTree::TElemSeq elems;
1231 ebbTree->getElementsInSphere( center, radius, elems );
1233 foundElems.insert( foundElems.end(), elems.begin(), elems.end() );
1236 //=======================================================================
1238 * Return elements whose bounding box intersects a given bounding box
1240 //=======================================================================
1242 void SMESH_ElementSearcherImpl::
1243 GetElementsInBox( const Bnd_B3d& box,
1244 SMDSAbs_ElementType type,
1245 std::vector< const SMDS_MeshElement* >& foundElems)
1247 _elementType = type;
1248 ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
1250 ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt, getTolerance() );
1252 ElementBndBoxTree::TElemSeq elems;
1253 ebbTree->getElementsInBox( box, elems );
1255 foundElems.insert( foundElems.end(), elems.begin(), elems.end() );
1258 //=======================================================================
1260 * \brief Return a projection of a given point to a mesh.
1261 * Optionally return the closest element
1263 //=======================================================================
1265 gp_XYZ SMESH_ElementSearcherImpl::Project(const gp_Pnt& point,
1266 SMDSAbs_ElementType type,
1267 const SMDS_MeshElement** closestElem)
1269 _elementType = type;
1270 if ( _mesh->GetMeshInfo().NbElements( _elementType ) == 0 )
1271 throw SALOME_Exception( LOCALIZED( "No elements of given type in the mesh" ));
1273 ElementBndBoxTree*& ebbTree = _ebbTree[ _elementType ];
1275 ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
1277 gp_XYZ p = point.XYZ();
1278 ElementBndBoxTree* ebbLeaf = ebbTree->getLeafAtPoint( p );
1279 const Bnd_B3d* box = ebbLeaf ? ebbLeaf->getBox() : ebbTree->getBox();
1280 gp_XYZ pMin = box->CornerMin(), pMax = box->CornerMax();
1281 double radius = Precision::Infinite();
1282 if ( ebbLeaf || !box->IsOut( p ))
1284 for ( int i = 1; i <= 3; ++i )
1286 double d = 0.5 * ( pMax.Coord(i) - pMin.Coord(i) );
1287 if ( d > Precision::Confusion() )
1288 radius = Min( d, radius );
1291 radius /= ebbTree->getHeight( /*full=*/true );
1293 else // p outside of box
1295 for ( int i = 1; i <= 3; ++i )
1298 if ( point.Coord(i) < pMin.Coord(i) )
1299 d = pMin.Coord(i) - point.Coord(i);
1300 else if ( point.Coord(i) > pMax.Coord(i) )
1301 d = point.Coord(i) - pMax.Coord(i);
1302 if ( d > Precision::Confusion() )
1303 radius = Min( d, radius );
1307 ElementBndBoxTree::TElemSeq elems;
1308 ebbTree->getElementsInSphere( p, radius, elems );
1309 while ( elems.empty() && radius < 1e100 )
1312 ebbTree->getElementsInSphere( p, radius, elems );
1314 gp_XYZ proj, bestProj;
1315 const SMDS_MeshElement* elem = 0;
1316 double minDist = Precision::Infinite();
1317 ElementBndBoxTree::TElemSeq::iterator e = elems.begin();
1318 for ( ; e != elems.end(); ++e )
1320 double d = SMESH_MeshAlgos::GetDistance( *e, point, &proj );
1328 if ( minDist > radius )
1330 ElementBndBoxTree::TElemSeq elems2;
1331 ebbTree->getElementsInSphere( p, minDist, elems2 );
1332 for ( e = elems2.begin(); e != elems2.end(); ++e )
1334 if ( elems.count( *e ))
1336 double d = SMESH_MeshAlgos::GetDistance( *e, point, &proj );
1345 if ( closestElem ) *closestElem = elem;
1350 //=======================================================================
1352 * \brief Return true if the point is IN or ON of the element
1354 //=======================================================================
1356 bool SMESH_MeshAlgos::IsOut( const SMDS_MeshElement* element, const gp_Pnt& point, double tol )
1358 if ( element->GetType() == SMDSAbs_Volume)
1360 return SMDS_VolumeTool( element ).IsOut( point.X(), point.Y(), point.Z(), tol );
1363 // get ordered nodes
1365 std::vector< SMESH_TNodeXYZ > xyz; xyz.reserve( element->NbNodes()+1 );
1367 SMDS_NodeIteratorPtr nodeIt = element->interlacedNodesIterator();
1368 for ( int i = 0; nodeIt->more(); ++i )
1369 xyz.push_back( SMESH_TNodeXYZ( nodeIt->next() ));
1371 int i, nbNodes = (int) xyz.size(); // central node of biquadratic is missing
1373 if ( element->GetType() == SMDSAbs_Face ) // --------------------------------------------------
1375 // compute face normal
1376 gp_Vec faceNorm(0,0,0);
1377 xyz.push_back( xyz.front() );
1378 for ( i = 0; i < nbNodes; ++i )
1380 gp_Vec edge1( xyz[i+1], xyz[i]);
1381 gp_Vec edge2( xyz[i+1], xyz[(i+2)%nbNodes] );
1382 faceNorm += edge1 ^ edge2;
1384 double fNormSize = faceNorm.Magnitude();
1385 if ( fNormSize <= tol )
1387 // degenerated face: point is out if it is out of all face edges
1388 for ( i = 0; i < nbNodes; ++i )
1390 SMDS_LinearEdge edge( xyz[i]._node, xyz[i+1]._node );
1391 if ( !IsOut( &edge, point, tol ))
1396 faceNorm /= fNormSize;
1398 // check if the point lays on face plane
1399 gp_Vec n2p( xyz[0], point );
1400 double dot = n2p * faceNorm;
1401 if ( Abs( dot ) > tol ) // not on face plane
1404 if ( nbNodes > 3 ) // maybe the face is not planar
1406 double elemThick = 0;
1407 for ( i = 1; i < nbNodes; ++i )
1409 gp_Vec n2n( xyz[0], xyz[i] );
1410 elemThick = Max( elemThick, Abs( n2n * faceNorm ));
1412 isOut = Abs( dot ) > elemThick + tol;
1418 // check if point is out of face boundary:
1419 // define it by closest transition of a ray point->infinity through face boundary
1420 // on the face plane.
1421 // First, find normal of a plane perpendicular to face plane, to be used as a cutting tool
1422 // to find intersections of the ray with the boundary.
1424 gp_Vec plnNorm = ray ^ faceNorm;
1425 double n2pSize = plnNorm.Magnitude();
1426 if ( n2pSize <= tol ) return false; // point coincides with the first node
1427 if ( n2pSize * n2pSize > fNormSize * 100 ) return true; // point is very far
1429 // for each node of the face, compute its signed distance to the cutting plane
1430 std::vector<double> dist( nbNodes + 1);
1431 for ( i = 0; i < nbNodes; ++i )
1433 gp_Vec n2p( xyz[i], point );
1434 dist[i] = n2p * plnNorm;
1436 dist.back() = dist.front();
1437 // find the closest intersection
1439 double rClosest = 0, distClosest = 1e100;
1441 for ( i = 0; i < nbNodes; ++i )
1444 if ( fabs( dist[i] ) < tol )
1446 else if ( fabs( dist[i+1]) < tol )
1448 else if ( dist[i] * dist[i+1] < 0 )
1449 r = dist[i] / ( dist[i] - dist[i+1] );
1451 continue; // no intersection
1452 gp_Pnt pInt = xyz[i] * (1.-r) + xyz[i+1] * r;
1453 gp_Vec p2int( point, pInt);
1454 double intDist = p2int.SquareMagnitude();
1455 if ( intDist < distClosest )
1460 distClosest = intDist;
1464 return true; // no intesections - out
1466 // analyse transition
1467 gp_Vec edge( xyz[iClosest], xyz[iClosest+1] );
1468 gp_Vec edgeNorm = -( edge ^ faceNorm ); // normal to intersected edge pointing out of face
1469 gp_Vec p2int ( point, pClosest );
1470 bool out = (edgeNorm * p2int) < -tol;
1471 if ( rClosest > 0. && rClosest < 1. ) // not node intersection
1474 // the ray passes through a face node; analyze transition through an adjacent edge
1475 gp_Pnt p1 = xyz[ (rClosest == 0.) ? ((iClosest+nbNodes-1) % nbNodes) : (iClosest+1) ];
1476 gp_Pnt p2 = xyz[ (rClosest == 0.) ? iClosest : ((iClosest+2) % nbNodes) ];
1477 gp_Vec edgeAdjacent( p1, p2 );
1478 gp_Vec edgeNorm2 = -( edgeAdjacent ^ faceNorm );
1479 bool out2 = (edgeNorm2 * p2int) < -tol;
1481 bool covexCorner = ( edgeNorm * edgeAdjacent * (rClosest==1. ? 1. : -1.)) < 0;
1482 return covexCorner ? (out || out2) : (out && out2);
1485 if ( element->GetType() == SMDSAbs_Edge ) // --------------------------------------------------
1487 // point is out of edge if it is NOT ON any straight part of edge
1488 // (we consider quadratic edge as being composed of two straight parts)
1489 for ( i = 1; i < nbNodes; ++i )
1491 gp_Vec edge( xyz[i-1], xyz[i] );
1492 gp_Vec n1p ( xyz[i-1], point );
1493 double u = ( edge * n1p ) / edge.SquareMagnitude(); // param [0,1] on the edge
1495 if ( n1p.SquareMagnitude() < tol * tol )
1500 if ( point.SquareDistance( xyz[i] ) < tol * tol )
1504 gp_XYZ proj = ( 1. - u ) * xyz[i-1] + u * xyz[i]; // projection of the point on the edge
1505 double dist2 = point.SquareDistance( proj );
1506 if ( dist2 > tol * tol )
1508 return false; // point is ON this part
1513 // Node or 0D element -------------------------------------------------------------------------
1515 gp_Vec n2p ( xyz[0], point );
1516 return n2p.SquareMagnitude() > tol * tol;
1521 //=======================================================================
1524 // Position of a point relative to a segment
1528 // VERTEX 1 o----ON-----> VERTEX 2
1532 enum PositionName { POS_LEFT = 1, POS_VERTEX = 2, POS_RIGHT = 4, //POS_ON = 8,
1533 POS_ALL = POS_LEFT | POS_RIGHT | POS_VERTEX,
1534 POS_MAX = POS_RIGHT };
1538 int _index; // index of vertex or segment
1540 PointPos( PositionName n, int i=-1 ): _name(n), _index(i) {}
1541 bool operator < (const PointPos& other ) const
1543 if ( _name == other._name )
1544 return ( _index < 0 || other._index < 0 ) ? false : _index < other._index;
1545 return _name < other._name;
1549 //================================================================================
1551 * \brief Return position of a point relative to a segment
1552 * \param point2D - the point to analyze position of
1553 * \param segEnds - end points of segments
1554 * \param index0 - 0-based index of the first point of segment
1555 * \param posToFindOut - flags of positions to detect
1556 * \retval PointPos - point position
1558 //================================================================================
1560 PointPos getPointPosition( const gp_XY& point2D,
1561 const gp_XY* segEnds,
1562 const int index0 = 0,
1563 const int posToFindOut = POS_ALL)
1565 const gp_XY& p1 = segEnds[ index0 ];
1566 const gp_XY& p2 = segEnds[ index0+1 ];
1567 const gp_XY grad = p2 - p1;
1569 if ( posToFindOut & POS_VERTEX )
1571 // check if the point2D is at "vertex 1" zone
1572 gp_XY pp1[2] = { p1, gp_XY( p1.X() - grad.Y(),
1573 p1.Y() + grad.X() ) };
1574 if ( getPointPosition( point2D, pp1, 0, POS_LEFT|POS_RIGHT )._name == POS_LEFT )
1575 return PointPos( POS_VERTEX, index0 );
1577 // check if the point2D is at "vertex 2" zone
1578 gp_XY pp2[2] = { p2, gp_XY( p2.X() - grad.Y(),
1579 p2.Y() + grad.X() ) };
1580 if ( getPointPosition( point2D, pp2, 0, POS_LEFT|POS_RIGHT )._name == POS_RIGHT )
1581 return PointPos( POS_VERTEX, index0 + 1);
1583 double edgeEquation =
1584 ( point2D.X() - p1.X() ) * grad.Y() - ( point2D.Y() - p1.Y() ) * grad.X();
1585 return PointPos( edgeEquation < 0 ? POS_LEFT : POS_RIGHT, index0 );
1589 //=======================================================================
1591 * \brief Return minimal distance from a point to an element
1593 * Currently we ignore non-planarity and 2nd order of face
1595 //=======================================================================
1597 double SMESH_MeshAlgos::GetDistance( const SMDS_MeshElement* elem,
1598 const gp_Pnt& point,
1599 gp_XYZ* closestPnt )
1601 switch ( elem->GetType() )
1603 case SMDSAbs_Volume:
1604 return GetDistance( static_cast<const SMDS_MeshVolume*>( elem ), point, closestPnt );
1606 return GetDistance( static_cast<const SMDS_MeshFace*>( elem ), point, closestPnt );
1608 return GetDistance( static_cast<const SMDS_MeshEdge*>( elem ), point, closestPnt );
1610 if ( closestPnt ) *closestPnt = SMESH_TNodeXYZ( elem );
1611 return point.Distance( SMESH_TNodeXYZ( elem ));
1617 //=======================================================================
1619 * \brief Return minimal distance from a point to a face
1621 * Currently we ignore non-planarity and 2nd order of face
1623 //=======================================================================
1625 double SMESH_MeshAlgos::GetDistance( const SMDS_MeshFace* face,
1626 const gp_Pnt& point,
1627 gp_XYZ* closestPnt )
1629 const double badDistance = -1;
1630 if ( !face ) return badDistance;
1632 int nbCorners = face->NbCornerNodes();
1633 if ( nbCorners > 3 )
1635 std::vector< const SMDS_MeshNode* > nodes;
1636 int nbTria = SMESH_MeshAlgos::Triangulate().GetTriangles( face, nodes );
1638 double minDist = Precision::Infinite();
1640 for ( int i = 0; i < 3 * nbTria; i += 3 )
1642 SMDS_FaceOfNodes triangle( nodes[i], nodes[i+1], nodes[i+2] );
1643 double dist = GetDistance( &triangle, point, closestPnt );
1644 if ( dist < minDist )
1657 // coordinates of nodes (medium nodes, if any, ignored)
1658 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
1659 std::vector<gp_XYZ> xyz( TXyzIterator( face->nodesIterator()), TXyzIterator() );
1662 // transformation to get xyz[0] lies on the origin, xyz[1] lies on the Z axis,
1663 // and xyz[2] lies in the XZ plane. This is to pass to 2D space on XZ plane.
1665 gp_Vec OZ ( xyz[0], xyz[1] );
1666 gp_Vec OX ( xyz[0], xyz[2] );
1667 if ( OZ.Magnitude() < std::numeric_limits<double>::min() )
1669 if ( xyz.size() < 4 ) return badDistance;
1670 OZ = gp_Vec ( xyz[0], xyz[2] );
1671 OX = gp_Vec ( xyz[0], xyz[3] );
1675 tgtCS = gp_Ax3( xyz[0], OZ, OX );
1677 catch ( Standard_Failure& ) {
1680 trsf.SetTransformation( tgtCS );
1682 // move all the nodes to 2D
1683 std::vector<gp_XY> xy( xyz.size() );
1684 for ( size_t i = 0; i < 3; ++i )
1686 gp_XYZ p3d = xyz[i];
1687 trsf.Transforms( p3d );
1688 xy[i].SetCoord( p3d.X(), p3d.Z() );
1690 xyz.back() = xyz.front();
1691 xy.back() = xy.front();
1693 // // move the point in 2D
1694 gp_XYZ tmpPnt = point.XYZ();
1695 trsf.Transforms( tmpPnt );
1696 gp_XY point2D( tmpPnt.X(), tmpPnt.Z() );
1698 // loop on edges of the face to analyze point position ralative to the face
1699 std::vector< PointPos > pntPosByType[ POS_MAX + 1 ];
1700 for ( size_t i = 1; i < xy.size(); ++i )
1702 PointPos pos = getPointPosition( point2D, &xy[0], i-1 );
1703 pntPosByType[ pos._name ].push_back( pos );
1708 double dist = badDistance;
1710 if ( pntPosByType[ POS_LEFT ].size() > 0 ) // point is most close to an edge
1712 PointPos& pos = pntPosByType[ POS_LEFT ][0];
1714 gp_Vec edge( xyz[ pos._index ], xyz[ pos._index+1 ]);
1715 gp_Vec n1p ( xyz[ pos._index ], point );
1716 double u = ( edge * n1p ) / edge.SquareMagnitude(); // param [0,1] on the edge
1717 gp_XYZ proj = xyz[ pos._index ] + u * edge.XYZ(); // projection on the edge
1718 dist = point.Distance( proj );
1719 if ( closestPnt ) *closestPnt = proj;
1722 else if ( pntPosByType[ POS_RIGHT ].size() >= 2 ) // point is inside the face
1724 dist = Abs( tmpPnt.Y() );
1727 if ( dist < std::numeric_limits<double>::min() ) {
1728 *closestPnt = point.XYZ();
1732 trsf.Inverted().Transforms( tmpPnt );
1733 *closestPnt = tmpPnt;
1738 else if ( pntPosByType[ POS_VERTEX ].size() > 0 ) // point is most close to a node
1740 double minDist2 = Precision::Infinite();
1741 for ( size_t i = 0; i < pntPosByType[ POS_VERTEX ].size(); ++i )
1743 PointPos& pos = pntPosByType[ POS_VERTEX ][i];
1745 double d2 = point.SquareDistance( xyz[ pos._index ]);
1746 if ( minDist2 > d2 )
1748 if ( closestPnt ) *closestPnt = xyz[ pos._index ];
1752 dist = Sqrt( minDist2 );
1758 //=======================================================================
1760 * \brief Return minimal distance from a point to an edge
1762 //=======================================================================
1764 double SMESH_MeshAlgos::GetDistance( const SMDS_MeshEdge* seg,
1765 const gp_Pnt& point,
1766 gp_XYZ* closestPnt )
1768 double dist = Precision::Infinite();
1769 if ( !seg ) return dist;
1771 int i = 0, nbNodes = seg->NbNodes();
1773 std::vector< SMESH_TNodeXYZ > xyz( nbNodes );
1774 for ( SMDS_NodeIteratorPtr nodeIt = seg->interlacedNodesIterator(); nodeIt->more(); i++ )
1775 xyz[ i ].Set( nodeIt->next() );
1777 for ( i = 1; i < nbNodes; ++i )
1779 gp_Vec edge( xyz[i-1], xyz[i] );
1780 gp_Vec n1p ( xyz[i-1], point );
1781 double d, u = ( edge * n1p ) / edge.SquareMagnitude(); // param [0,1] on the edge
1783 if (( d = n1p.SquareMagnitude() ) < dist ) {
1785 if ( closestPnt ) *closestPnt = xyz[i-1];
1788 else if ( u >= 1. ) {
1789 if (( d = point.SquareDistance( xyz[i] )) < dist ) {
1791 if ( closestPnt ) *closestPnt = xyz[i];
1795 gp_XYZ proj = xyz[i-1] + u * edge.XYZ(); // projection of the point on the edge
1796 if (( d = point.SquareDistance( proj )) < dist ) {
1798 if ( closestPnt ) *closestPnt = proj;
1802 return Sqrt( dist );
1805 //=======================================================================
1807 * \brief Return minimal distance from a point to a volume
1809 * Currently we ignore non-planarity and 2nd order
1811 //=======================================================================
1813 double SMESH_MeshAlgos::GetDistance( const SMDS_MeshVolume* volume,
1814 const gp_Pnt& point,
1815 gp_XYZ* closestPnt )
1817 SMDS_VolumeTool vTool( volume );
1818 vTool.SetExternalNormal();
1819 const int iQ = volume->IsQuadratic() ? 2 : 1;
1822 double minDist = 1e100, dist;
1823 gp_XYZ closeP = point.XYZ();
1825 for ( int iF = 0; iF < vTool.NbFaces(); ++iF )
1827 // skip a facet with normal not "looking at" the point
1828 if ( !vTool.GetFaceNormal( iF, n[0], n[1], n[2] ) ||
1829 !vTool.GetFaceBaryCenter( iF, bc[0], bc[1], bc[2] ))
1831 gp_XYZ bcp = point.XYZ() - gp_XYZ( bc[0], bc[1], bc[2] );
1832 if ( gp_XYZ( n[0], n[1], n[2] ) * bcp < -1e-12 )
1835 // find distance to a facet
1836 const SMDS_MeshNode** nodes = vTool.GetFaceNodes( iF );
1837 switch ( vTool.NbFaceNodes( iF ) / iQ ) {
1840 SMDS_FaceOfNodes tmpFace( nodes[0], nodes[ 1*iQ ], nodes[ 2*iQ ] );
1841 dist = GetDistance( &tmpFace, point, closestPnt );
1846 SMDS_FaceOfNodes tmpFace( nodes[0], nodes[ 1*iQ ], nodes[ 2*iQ ], nodes[ 3*iQ ]);
1847 dist = GetDistance( &tmpFace, point, closestPnt );
1851 std::vector<const SMDS_MeshNode *> nvec( nodes, nodes + vTool.NbFaceNodes( iF ));
1852 SMDS_PolygonalFaceOfNodes tmpFace( nvec );
1853 dist = GetDistance( &tmpFace, point, closestPnt );
1855 if ( dist < minDist )
1859 if ( closestPnt ) closeP = *closestPnt;
1864 if ( closestPnt ) *closestPnt = closeP;
1868 return 0; // point is inside the volume
1871 //================================================================================
1873 * \brief Returns barycentric coordinates of a point within a triangle.
1874 * A not returned bc2 = 1. - bc0 - bc1.
1875 * The point lies within the triangle if ( bc0 >= 0 && bc1 >= 0 && bc0+bc1 <= 1 )
1877 //================================================================================
1879 void SMESH_MeshAlgos::GetBarycentricCoords( const gp_XY& p,
1886 const double // matrix 2x2
1887 T11 = t0.X()-t2.X(), T12 = t1.X()-t2.X(),
1888 T21 = t0.Y()-t2.Y(), T22 = t1.Y()-t2.Y();
1889 const double Tdet = T11*T22 - T12*T21; // matrix determinant
1890 if ( Abs( Tdet ) < std::numeric_limits<double>::min() )
1896 const double t11 = T22, t12 = -T12, t21 = -T21, t22 = T11;
1898 const double r11 = p.X()-t2.X(), r12 = p.Y()-t2.Y();
1899 // barycentric coordinates: multiply matrix by vector
1900 bc0 = (t11 * r11 + t12 * r12)/Tdet;
1901 bc1 = (t21 * r11 + t22 * r12)/Tdet;
1904 //=======================================================================
1905 //function : FindFaceInSet
1906 //purpose : Return a face having linked nodes n1 and n2 and which is
1907 // - not in avoidSet,
1908 // - in elemSet provided that !elemSet.empty()
1909 // i1 and i2 optionally returns indices of n1 and n2
1910 //=======================================================================
1912 const SMDS_MeshElement*
1913 SMESH_MeshAlgos::FindFaceInSet(const SMDS_MeshNode* n1,
1914 const SMDS_MeshNode* n2,
1915 const TIDSortedElemSet& elemSet,
1916 const TIDSortedElemSet& avoidSet,
1922 const SMDS_MeshElement* face = 0;
1924 SMDS_ElemIteratorPtr invElemIt = n1->GetInverseElementIterator(SMDSAbs_Face);
1925 while ( invElemIt->more() && !face ) // loop on inverse faces of n1
1927 const SMDS_MeshElement* elem = invElemIt->next();
1928 if (avoidSet.count( elem ))
1930 if ( !elemSet.empty() && !elemSet.count( elem ))
1933 i1 = elem->GetNodeIndex( n1 );
1934 // find a n2 linked to n1
1935 int nbN = elem->IsQuadratic() ? elem->NbNodes()/2 : elem->NbNodes();
1936 for ( int di = -1; di < 2 && !face; di += 2 )
1938 i2 = (i1+di+nbN) % nbN;
1939 if ( elem->GetNode( i2 ) == n2 )
1942 if ( !face && elem->IsQuadratic())
1944 // analysis for quadratic elements using all nodes
1945 SMDS_NodeIteratorPtr anIter = elem->interlacedNodesIterator();
1946 const SMDS_MeshNode* prevN = static_cast<const SMDS_MeshNode*>( anIter->next() );
1947 for ( i1 = -1, i2 = 0; anIter->more() && !face; i1++, i2++ )
1949 const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( anIter->next() );
1950 if ( n1 == prevN && n2 == n )
1954 else if ( n2 == prevN && n1 == n )
1956 face = elem; std::swap( i1, i2 );
1962 if ( n1ind ) *n1ind = i1;
1963 if ( n2ind ) *n2ind = i2;
1967 //================================================================================
1969 * Return sharp edges of faces and non-manifold ones. Optionally adds existing edges.
1971 //================================================================================
1973 std::vector< SMESH_MeshAlgos::Edge >
1974 SMESH_MeshAlgos::FindSharpEdges( SMDS_Mesh* theMesh,
1976 bool theAddExisting )
1978 std::vector< Edge > resultEdges;
1979 if ( !theMesh ) return resultEdges;
1981 typedef std::pair< bool, const SMDS_MeshNode* > TIsSharpAndMedium;
1982 typedef NCollection_DataMap< SMESH_TLink, TIsSharpAndMedium, SMESH_TLink > TLinkSharpMap;
1984 TLinkSharpMap linkIsSharp( theMesh->NbFaces() );
1985 TIsSharpAndMedium sharpMedium( true, 0 );
1986 bool & isSharp = sharpMedium.first;
1987 const SMDS_MeshNode* & nMedium = sharpMedium.second;
1989 if ( theAddExisting )
1991 for ( SMDS_EdgeIteratorPtr edgeIt = theMesh->edgesIterator(); edgeIt->more(); )
1993 const SMDS_MeshElement* edge = edgeIt->next();
1994 nMedium = ( edge->IsQuadratic() ) ? edge->GetNode(2) : 0;
1995 linkIsSharp.Bind( SMESH_TLink( edge->GetNode(0), edge->GetNode(1)), sharpMedium );
1999 // check angles between face normals
2001 const double angleCos = Cos( theAngle * M_PI / 180. ), angleCos2 = angleCos * angleCos;
2002 gp_XYZ norm1, norm2;
2003 std::vector< const SMDS_MeshNode* > faceNodes, linkNodes(2);
2004 std::vector<const SMDS_MeshElement *> linkFaces;
2006 int nbSharp = linkIsSharp.Extent();
2007 for ( SMDS_FaceIteratorPtr faceIt = theMesh->facesIterator(); faceIt->more(); )
2009 const SMDS_MeshElement* face = faceIt->next();
2010 size_t nbCorners = face->NbCornerNodes();
2012 faceNodes.assign( face->begin_nodes(), face->end_nodes() );
2013 if ( faceNodes.size() == nbCorners )
2014 faceNodes.resize( nbCorners * 2, 0 );
2016 const SMDS_MeshNode* nPrev = faceNodes[ nbCorners-1 ];
2017 for ( size_t i = 0; i < nbCorners; ++i )
2019 SMESH_TLink link( nPrev, faceNodes[i] );
2020 if ( !linkIsSharp.IsBound( link ))
2022 linkNodes[0] = link.node1();
2023 linkNodes[1] = link.node2();
2025 theMesh->GetElementsByNodes( linkNodes, linkFaces, SMDSAbs_Face );
2028 if ( linkFaces.size() > 2 )
2032 else if ( linkFaces.size() == 2 &&
2033 FaceNormal( linkFaces[0], norm1, /*normalize=*/false ) &&
2034 FaceNormal( linkFaces[1], norm2, /*normalize=*/false ))
2036 double dot = norm1 * norm2; // == cos * |norm1| * |norm2|
2037 if (( dot < 0 ) == ( angleCos < 0 ))
2039 double cos2 = dot * dot / norm1.SquareModulus() / norm2.SquareModulus();
2040 isSharp = ( angleCos < 0 ) ? ( cos2 > angleCos2 ) : ( cos2 < angleCos2 );
2044 isSharp = ( angleCos > 0 );
2047 nMedium = faceNodes[( i-1+nbCorners ) % nbCorners + nbCorners ];
2049 linkIsSharp.Bind( link, sharpMedium );
2053 nPrev = faceNodes[i];
2057 resultEdges.resize( nbSharp );
2058 TLinkSharpMap::Iterator linkIsSharpIter( linkIsSharp );
2059 for ( int i = 0; linkIsSharpIter.More() && i < nbSharp; linkIsSharpIter.Next() )
2061 const SMESH_TLink& link = linkIsSharpIter.Key();
2062 const TIsSharpAndMedium& isSharpMedium = linkIsSharpIter.Value();
2063 if ( isSharpMedium.first )
2065 Edge & edge = resultEdges[ i++ ];
2066 edge._node1 = link.node1();
2067 edge._node2 = link.node2();
2068 edge._medium = isSharpMedium.second;
2075 //================================================================================
2077 * Distribute all faces of the mesh between groups using given edges as group boundaries
2079 //================================================================================
2081 std::vector< std::vector< const SMDS_MeshElement* > >
2082 SMESH_MeshAlgos::SeparateFacesByEdges( SMDS_Mesh* theMesh, const std::vector< Edge >& theEdges )
2084 std::vector< std::vector< const SMDS_MeshElement* > > groups;
2085 if ( !theMesh ) return groups;
2087 // build map of face edges (SMESH_TLink) and their faces
2089 typedef std::vector< const SMDS_MeshElement* > TFaceVec;
2090 typedef NCollection_DataMap< SMESH_TLink, TFaceVec, SMESH_TLink > TFacesByLinks;
2091 TFacesByLinks facesByLink( theMesh->NbFaces() );
2093 std::vector< const SMDS_MeshNode* > faceNodes;
2094 for ( SMDS_FaceIteratorPtr faceIt = theMesh->facesIterator(); faceIt->more(); )
2096 const SMDS_MeshElement* face = faceIt->next();
2097 size_t nbCorners = face->NbCornerNodes();
2099 faceNodes.assign( face->begin_nodes(), face->end_nodes() );
2100 faceNodes.resize( nbCorners + 1 );
2101 faceNodes[ nbCorners ] = faceNodes[0];
2103 face->setIsMarked( false );
2105 for ( size_t i = 0; i < nbCorners; ++i )
2107 SMESH_TLink link( faceNodes[i], faceNodes[i+1] );
2108 TFaceVec* linkFaces = facesByLink.ChangeSeek( link );
2111 linkFaces = facesByLink.Bound( link, TFaceVec() );
2112 linkFaces->reserve(2);
2114 linkFaces->push_back( face );
2118 // remove the given edges from facesByLink map
2120 for ( size_t i = 0; i < theEdges.size(); ++i )
2122 SMESH_TLink link( theEdges[i]._node1, theEdges[i]._node2 );
2123 facesByLink.UnBind( link );
2126 // faces connected via links of facesByLink map form a group
2128 while ( !facesByLink.IsEmpty() )
2130 groups.push_back( TFaceVec() );
2131 TFaceVec & group = groups.back();
2133 group.push_back( TFacesByLinks::Iterator( facesByLink ).Value()[0] );
2134 group.back()->setIsMarked( true );
2136 for ( size_t iF = 0; iF < group.size(); ++iF )
2138 const SMDS_MeshElement* face = group[iF];
2139 size_t nbCorners = face->NbCornerNodes();
2140 faceNodes.assign( face->begin_nodes(), face->end_nodes() );
2141 faceNodes.resize( nbCorners + 1 );
2142 faceNodes[ nbCorners ] = faceNodes[0];
2144 for ( size_t iN = 0; iN < nbCorners; ++iN )
2146 SMESH_TLink link( faceNodes[iN], faceNodes[iN+1] );
2147 if ( const TFaceVec* faces = facesByLink.Seek( link ))
2149 const TFaceVec& faceNeighbors = *faces;
2150 for ( size_t i = 0; i < faceNeighbors.size(); ++i )
2151 if ( !faceNeighbors[i]->isMarked() )
2153 group.push_back( faceNeighbors[i] );
2154 faceNeighbors[i]->setIsMarked( true );
2156 facesByLink.UnBind( link );
2162 // find faces that are alone in its group; they were not in facesByLink
2165 for ( size_t i = 0; i < groups.size(); ++i )
2166 nbInGroups += groups[i].size();
2167 if ( nbInGroups < theMesh->NbFaces() )
2169 for ( SMDS_FaceIteratorPtr faceIt = theMesh->facesIterator(); faceIt->more(); )
2171 const SMDS_MeshElement* face = faceIt->next();
2172 if ( !face->isMarked() )
2174 groups.push_back( TFaceVec() );
2175 groups.back().push_back( face );
2183 //================================================================================
2185 * \brief Calculate normal of a mesh face
2187 //================================================================================
2189 bool SMESH_MeshAlgos::FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool normalized)
2191 if ( !F || F->GetType() != SMDSAbs_Face )
2194 normal.SetCoord(0,0,0);
2195 int nbNodes = F->NbCornerNodes();
2196 for ( int i = 0; i < nbNodes-2; ++i )
2199 for ( int n = 0; n < 3; ++n )
2201 const SMDS_MeshNode* node = F->GetNode( i + n );
2202 p[n].SetCoord( node->X(), node->Y(), node->Z() );
2204 normal += ( p[2] - p[1] ) ^ ( p[0] - p[1] );
2206 double size2 = normal.SquareModulus();
2207 bool ok = ( size2 > std::numeric_limits<double>::min() * std::numeric_limits<double>::min());
2208 if ( normalized && ok )
2209 normal /= sqrt( size2 );
2214 //================================================================================
2216 * \brief Return nodes common to two elements
2218 //================================================================================
2220 int SMESH_MeshAlgos::NbCommonNodes(const SMDS_MeshElement* e1,
2221 const SMDS_MeshElement* e2)
2224 for ( int i = 0 ; i < e1->NbNodes(); ++i )
2225 nb += ( e2->GetNodeIndex( e1->GetNode( i )) >= 0 );
2229 //================================================================================
2231 * \brief Return nodes common to two elements
2233 //================================================================================
2235 std::vector< const SMDS_MeshNode*> SMESH_MeshAlgos::GetCommonNodes(const SMDS_MeshElement* e1,
2236 const SMDS_MeshElement* e2)
2238 std::vector< const SMDS_MeshNode*> common;
2239 for ( int i = 0 ; i < e1->NbNodes(); ++i )
2240 if ( e2->GetNodeIndex( e1->GetNode( i )) >= 0 )
2241 common.push_back( e1->GetNode( i ));
2245 //================================================================================
2247 * \brief Return true if node1 encounters first in the face and node2, after
2249 //================================================================================
2251 bool SMESH_MeshAlgos::IsRightOrder( const SMDS_MeshElement* face,
2252 const SMDS_MeshNode* node0,
2253 const SMDS_MeshNode* node1 )
2255 int i0 = face->GetNodeIndex( node0 );
2256 int i1 = face->GetNodeIndex( node1 );
2257 if ( face->IsQuadratic() )
2259 if ( face->IsMediumNode( node0 ))
2261 i0 -= ( face->NbNodes()/2 - 1 );
2266 i1 -= ( face->NbNodes()/2 - 1 );
2271 return ( diff == 1 ) || ( diff == -face->NbNodes()+1 );
2274 //=======================================================================
2276 * \brief Partition given 1D elements into groups of contiguous edges.
2277 * A node where number of meeting edges != 2 is a group end.
2278 * An optional startNode is used to orient groups it belongs to.
2279 * \return a list of edge groups and a list of corresponding node groups.
2280 * If a group is closed, the first and last nodes of the group are same.
2282 //=======================================================================
2284 void SMESH_MeshAlgos::Get1DBranches( SMDS_ElemIteratorPtr theEdgeIt,
2285 TElemGroupVector& theEdgeGroups,
2286 TNodeGroupVector& theNodeGroups,
2287 const SMDS_MeshNode* theStartNode )
2292 // build map of nodes and their adjacent edges
2294 typedef std::vector< const SMDS_MeshNode* > TNodeVec;
2295 typedef std::vector< const SMDS_MeshElement* > TEdgeVec;
2296 typedef NCollection_DataMap< const SMDS_MeshNode*, TEdgeVec, SMESH_Hasher > TEdgesByNodeMap;
2297 TEdgesByNodeMap edgesByNode;
2299 while ( theEdgeIt->more() )
2301 const SMDS_MeshElement* edge = theEdgeIt->next();
2302 if ( edge->GetType() != SMDSAbs_Edge )
2305 const SMDS_MeshNode* nodes[2] = { edge->GetNode(0), edge->GetNode(1) };
2306 for ( int i = 0; i < 2; ++i )
2308 TEdgeVec* nodeEdges = edgesByNode.ChangeSeek( nodes[i] );
2311 nodeEdges = edgesByNode.Bound( nodes[i], TEdgeVec() );
2312 nodeEdges->reserve(2);
2314 nodeEdges->push_back( edge );
2318 if ( edgesByNode.IsEmpty() )
2322 // build edge branches
2324 TElemGroupVector branches(2);
2325 TNodeGroupVector nodeBranches(2);
2327 while ( !edgesByNode.IsEmpty() )
2329 if ( !theStartNode || !edgesByNode.IsBound( theStartNode ))
2331 theStartNode = TEdgesByNodeMap::Iterator( edgesByNode ).Key();
2334 size_t nbBranches = 0;
2335 bool startIsBranchEnd = false;
2337 while ( edgesByNode.IsBound( theStartNode ))
2339 // initialize a new branch
2342 if ( branches.size() < nbBranches )
2344 branches.push_back ( TEdgeVec() );
2345 nodeBranches.push_back( TNodeVec() );
2347 TEdgeVec & branch = branches [ nbBranches - 1 ];
2348 TNodeVec & nodeBranch = nodeBranches[ nbBranches - 1 ];
2352 TEdgeVec& edges = edgesByNode( theStartNode );
2353 startIsBranchEnd = ( edges.size() != 2 );
2356 const SMDS_MeshElement* startEdge = 0;
2357 for ( size_t i = 0; i < edges.size(); ++i )
2359 if ( !startEdge && edges[i] )
2361 startEdge = edges[i];
2364 nbEdges += bool( edges[i] );
2367 edgesByNode.UnBind( theStartNode );
2371 branch.push_back( startEdge );
2373 nodeBranch.push_back( theStartNode );
2374 nodeBranch.push_back( branch.back()->GetNode(0) );
2375 if ( nodeBranch.back() == theStartNode )
2376 nodeBranch.back() = branch.back()->GetNode(1);
2381 bool isBranchEnd = false;
2384 while (( !isBranchEnd ) && ( edgesPtr = edgesByNode.ChangeSeek( nodeBranch.back() )))
2386 TEdgeVec& edges = *edgesPtr;
2388 isBranchEnd = ( edges.size() != 2 );
2390 const SMDS_MeshNode* lastNode = nodeBranch.back();
2392 switch ( edges.size() )
2395 edgesByNode.UnBind( lastNode );
2400 if ( const SMDS_MeshElement* nextEdge = edges[ edges[0] == branch.back() ])
2402 branch.push_back( nextEdge );
2404 const SMDS_MeshNode* nextNode = nextEdge->GetNode(0);
2405 if ( nodeBranch.back() == nextNode )
2406 nextNode = nextEdge->GetNode(1);
2407 nodeBranch.push_back( nextNode );
2409 edgesByNode.UnBind( lastNode );
2415 for ( size_t i = 0; i < edges.size(); ++i )
2417 if ( edges[i] == branch.back() )
2419 nbEdges += bool( edges[i] );
2422 edgesByNode.UnBind( lastNode );
2425 } // while ( edgesByNode.IsBound( theStartNode ))
2428 // put the found branches to the result
2430 if ( nbBranches == 2 && !startIsBranchEnd ) // join two branches starting at the same node
2432 std::reverse( nodeBranches[0].begin(), nodeBranches[0].end() );
2433 nodeBranches[0].pop_back();
2434 nodeBranches[0].reserve( nodeBranches[0].size() + nodeBranches[1].size() );
2435 nodeBranches[0].insert( nodeBranches[0].end(),
2436 nodeBranches[1].begin(), nodeBranches[1].end() );
2438 std::reverse( branches[0].begin(), branches[0].end() );
2439 branches[0].reserve( branches[0].size() + branches[1].size() );
2440 branches[0].insert( branches[0].end(), branches[1].begin(), branches[1].end() );
2442 nodeBranches[1].clear();
2443 branches[1].clear();
2446 for ( size_t i = 0; i < nbBranches; ++i )
2448 if ( branches[i].empty() )
2451 theEdgeGroups.push_back( TEdgeVec() );
2452 theEdgeGroups.back().swap( branches[i] );
2454 theNodeGroups.push_back( TNodeVec() );
2455 theNodeGroups.back().swap( nodeBranches[i] );
2458 } // while ( !edgesByNode.IsEmpty() )
2463 //=======================================================================
2465 * \brief Return SMESH_NodeSearcher
2467 //=======================================================================
2469 SMESH_NodeSearcher* SMESH_MeshAlgos::GetNodeSearcher(SMDS_Mesh& mesh)
2471 return new SMESH_NodeSearcherImpl( &mesh );
2474 //=======================================================================
2476 * \brief Return SMESH_NodeSearcher
2478 //=======================================================================
2480 SMESH_NodeSearcher* SMESH_MeshAlgos::GetNodeSearcher(SMDS_ElemIteratorPtr elemIt)
2482 return new SMESH_NodeSearcherImpl( 0, elemIt );
2485 //=======================================================================
2487 * \brief Return SMESH_ElementSearcher
2489 //=======================================================================
2491 SMESH_ElementSearcher* SMESH_MeshAlgos::GetElementSearcher(SMDS_Mesh& mesh,
2494 return new SMESH_ElementSearcherImpl( mesh, tolerance );
2497 //=======================================================================
2499 * \brief Return SMESH_ElementSearcher acting on a sub-set of elements
2501 //=======================================================================
2503 SMESH_ElementSearcher* SMESH_MeshAlgos::GetElementSearcher(SMDS_Mesh& mesh,
2504 SMDS_ElemIteratorPtr elemIt,
2507 return new SMESH_ElementSearcherImpl( mesh, tolerance, elemIt );