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_Offset.cxx
23 // Created : Mon Dec 25 15:52:38 2017
24 // Author : Edward AGAPOV (eap)
26 #include "SMESH_MeshAlgos.hxx"
28 #include <SMDS_PolygonalFaceOfNodes.hxx>
29 #include "SMDS_Mesh.hxx"
31 #include <Utils_SALOME_Exception.hxx>
33 #include <Bnd_B3d.hxx>
34 #include <NCollection_Map.hxx>
38 #include <boost/container/flat_set.hpp>
39 #include <boost/dynamic_bitset.hpp>
43 const int theMaxNbFaces = 256; // max number of faces sharing a node
45 typedef NCollection_DataMap< const SMDS_MeshNode*, const SMDS_MeshNode*, SMESH_Hasher > TNNMap;
46 typedef NCollection_Map< SMESH_Link, SMESH_Link > TLinkMap;
48 //--------------------------------------------------------------------------------
50 * \brief Intersected face side storing a node created at this intersection
51 * and an intersected face
56 const SMDS_MeshNode* myNode[2]; // side nodes. WARNING: don't set them directly, use Set()
57 mutable SMESH_NodeXYZ myIntNode; // intersection node
58 const SMDS_MeshElement* myFace; // cutter face
59 int myIndex; // index of a node on the same link
61 CutLink(const SMDS_MeshNode* node1=0,
62 const SMDS_MeshNode* node2=0,
63 const SMDS_MeshElement* face=0,
64 const int index=0) { Set ( node1, node2, face, index ); }
66 void Set( const SMDS_MeshNode* node1,
67 const SMDS_MeshNode* node2,
68 const SMDS_MeshElement* face,
71 myNode[0] = node1; myNode[1] = node2; myFace = face; myIndex = index; myReverse = false;
72 if ( myNode[0] && ( myReverse = ( myNode[0]->GetID() > myNode[1]->GetID() )))
73 std::swap( myNode[0], myNode[1] );
75 const SMDS_MeshNode* IntNode() const { return myIntNode.Node(); }
76 const SMDS_MeshNode* Node1() const { return myNode[ myReverse ]; }
77 const SMDS_MeshNode* Node2() const { return myNode[ !myReverse ]; }
79 static Standard_Integer HashCode(const CutLink& link,
80 const Standard_Integer upper)
82 Standard_Integer n = ( link.myNode[0]->GetID() +
83 link.myNode[1]->GetID() +
85 return ::HashCode( n, upper );
87 static Standard_Boolean IsEqual(const CutLink& link1, const CutLink& link2 )
89 return ( link1.myNode[0] == link2.myNode[0] &&
90 link1.myNode[1] == link2.myNode[1] &&
91 link1.myIndex == link2.myIndex );
95 typedef NCollection_Map< CutLink, CutLink > TCutLinkMap;
97 //--------------------------------------------------------------------------------
99 * \brief Part of a divided face edge
103 const SMDS_MeshNode* myNode1;
104 const SMDS_MeshNode* myNode2;
105 int myIndex; // positive -> side index, negative -> State
106 const SMDS_MeshElement* myFace;
108 enum State { _INTERNAL = -1, _COPLANAR = -2, _PENDING = -3 };
110 void Set( const SMDS_MeshNode* Node1,
111 const SMDS_MeshNode* Node2,
112 const SMDS_MeshElement* Face = 0,
113 int EdgeIndex = _INTERNAL )
114 { myNode1 = Node1; myNode2 = Node2; myIndex = EdgeIndex; myFace = Face; }
116 // bool HasSameNode( const EdgePart& other ) { return ( myNode1 == other.myNode1 ||
117 // myNode1 == other.myNode2 ||
118 // myNode2 == other.myNode1 ||
119 // myNode2 == other.myNode2 );
121 bool IsInternal() const { return myIndex < 0; }
122 bool IsTwin( const EdgePart& e ) const { return myNode1 == e.myNode2 && myNode2 == e.myNode1; }
123 bool IsSame( const EdgePart& e ) const {
124 return (( myNode1 == e.myNode2 && myNode2 == e.myNode1 ) ||
125 ( myNode1 == e.myNode1 && myNode2 == e.myNode2 )); }
126 bool ReplaceCoplanar( const EdgePart& e );
127 operator SMESH_Link() const { return SMESH_Link( myNode1, myNode2 ); }
128 operator gp_Vec() const { return SMESH_NodeXYZ( myNode2 ) - SMESH_NodeXYZ( myNode1 ); }
131 //--------------------------------------------------------------------------------
133 * \brief Loop of EdgePart's forming a new face which is a part of CutFace
135 struct EdgeLoop : public SMDS_PolygonalFaceOfNodes
137 std::vector< const EdgePart* > myLinks;
138 bool myIsBndConnected; //!< is there a path to CutFace side edges
139 bool myHasPending; //!< an edge encounters twice
141 EdgeLoop() : SMDS_PolygonalFaceOfNodes( std::vector<const SMDS_MeshNode *>() ) {}
142 void Clear() { myLinks.clear(); myIsBndConnected = false; myHasPending = false; }
143 bool SetConnected() { bool was = myIsBndConnected; myIsBndConnected = true; return !was; }
144 size_t Contains( const SMDS_MeshNode* n ) const
146 for ( size_t i = 0; i < myLinks.size(); ++i )
147 if ( myLinks[i]->myNode1 == n ) return i + 1;
150 virtual int NbNodes() const { return myLinks.size(); }
151 virtual SMDS_ElemIteratorPtr nodesIterator() const
153 return setNodes(), SMDS_PolygonalFaceOfNodes::nodesIterator();
155 virtual SMDS_NodeIteratorPtr nodeIterator() const
157 return setNodes(), SMDS_PolygonalFaceOfNodes::nodeIterator();
159 void setNodes() const //!< set nodes to SMDS_PolygonalFaceOfNodes
161 EdgeLoop* me = const_cast<EdgeLoop*>( this );
162 me->myNodes.resize( NbNodes() );
164 for ( size_t i = 1; i < myNodes.size(); ++i ) {
165 if ( myLinks[ i ]->myNode1->GetID() < myLinks[ iMin ]->myNode1->GetID() )
168 for ( size_t i = 0; i < myNodes.size(); ++i )
169 me->myNodes[ i ] = myLinks[ ( iMin + i ) % myNodes.size() ]->myNode1;
173 //--------------------------------------------------------------------------------
175 * \brief Set of EdgeLoop's constructed from a CutFace
179 std::vector< EdgeLoop > myLoops; //!< buffer of EdgeLoop's
180 size_t myNbLoops; //!< number of constructed loops
182 const EdgePart* myEdge0; //!< & CutFace.myLinks[0]
183 size_t myNbUsedEdges; //!< nb of EdgePart's added to myLoops
184 boost::dynamic_bitset<> myIsUsedEdge; //!< is i-th EdgePart of CutFace is in any EdgeLoop
185 std::vector< EdgeLoop* > myLoopOfEdge; //!< EdgeLoop of CutFace.myLinks[i]
186 std::vector< EdgePart* > myCandidates; //!< EdgePart's starting at the same node
188 EdgeLoopSet(): myLoops(100) {}
190 void Init( const std::vector< EdgePart >& edges )
192 size_t nb = edges.size();
193 myEdge0 = & edges[0];
196 myIsUsedEdge.reset();
197 myIsUsedEdge.resize( nb, false );
198 myLoopOfEdge.clear();
199 myLoopOfEdge.resize( nb, (EdgeLoop*) 0 );
201 EdgeLoop& AddNewLoop()
203 if ( ++myNbLoops >= myLoops.size() )
204 myLoops.resize( myNbLoops + 10 );
205 myLoops[ myNbLoops-1 ].Clear();
206 return myLoops[ myNbLoops-1 ];
208 bool AllEdgesUsed() const { return myNbUsedEdges == myLoopOfEdge.size(); }
210 bool AddEdge( EdgePart& edge )
212 size_t i = Index( edge );
213 if ( myIsUsedEdge[ i ])
215 myLoops[ myNbLoops-1 ].myLinks.push_back( &edge );
216 myLoopOfEdge[ i ] = & myLoops[ myNbLoops-1 ];
217 myIsUsedEdge[ i ] = true;
221 void Erase( EdgeLoop* loop )
223 for ( size_t iE = 0; iE < loop->myLinks.size(); ++iE )
224 myLoopOfEdge[ Index( *loop->myLinks[ iE ] )] = 0;
227 void Join( EdgeLoop& loop1, size_t iAfterConcact,
228 EdgeLoop& loop2, size_t iFromEdge2 )
230 std::vector< const EdgePart* > linksAfterContact( loop1.myLinks.begin() + iAfterConcact,
231 loop1.myLinks.end() );
232 loop1.myLinks.reserve( loop2.myLinks.size() + loop1.myLinks.size() );
233 loop1.myLinks.resize( iAfterConcact );
234 loop1.myLinks.insert( loop1.myLinks.end(),
235 loop2.myLinks.begin() + iFromEdge2, loop2.myLinks.end() );
236 loop1.myLinks.insert( loop1.myLinks.end(),
237 loop2.myLinks.begin(), loop2.myLinks.begin() + iFromEdge2 );
238 loop1.myLinks.insert( loop1.myLinks.end(),
239 linksAfterContact.begin(), linksAfterContact.end() );
240 loop1.myIsBndConnected = loop2.myIsBndConnected;
242 for ( size_t iE = 0; iE < loop1.myLinks.size(); ++iE )
243 myLoopOfEdge[ Index( *loop1.myLinks[ iE ] )] = & loop1;
245 size_t Index( const EdgePart& edge ) const { return &edge - myEdge0; }
246 EdgeLoop* GetLoopOf( const EdgePart* edge ) { return myLoopOfEdge[ Index( *edge )]; }
249 //--------------------------------------------------------------------------------
251 * \brief Intersections of a face
255 mutable std::vector< EdgePart > myLinks;
256 const SMDS_MeshElement* myInitFace;
258 CutFace( const SMDS_MeshElement* face ): myInitFace( face ) {}
259 void AddEdge( const CutLink& p1,
261 const SMDS_MeshElement* cutter,
263 const int iNotOnPlane = -1) const;
264 void AddPoint( const CutLink& p1, const CutLink& p2, double tol ) const;
265 bool ReplaceNodes( const TNNMap& theRm2KeepMap ) const;
267 int NbInternalEdges() const;
268 void MakeLoops( EdgeLoopSet& loops, const gp_XYZ& theFaceNorm ) const;
269 bool RemoveInternalLoops( EdgeLoopSet& theLoops ) const;
270 void CutOffLoops( EdgeLoopSet& theLoops,
271 const double theSign,
272 const std::vector< gp_XYZ >& theNormals,
273 std::vector< EdgePart >& theCutOffLinks,
274 TLinkMap& theCutOffCoplanarLinks) const;
275 void InitLinks() const;
276 bool IsCoplanar( const EdgePart* edge ) const;
278 static Standard_Integer HashCode(const CutFace& f, const Standard_Integer upper)
280 return ::HashCode( f.myInitFace->GetID(), upper );
282 static Standard_Boolean IsEqual(const CutFace& f1, const CutFace& f2 )
284 return f1.myInitFace == f2.myInitFace;
290 EdgePart* getTwin( const EdgePart* edge ) const;
293 typedef NCollection_Map< CutFace, CutFace > TCutFaceMap;
295 //--------------------------------------------------------------------------------
297 * \brief Intersection point of two edges of co-planar triangles
301 size_t myEdgeInd[2]; //!< edge indices of triangles
302 double myU [2]; //!< parameter [0,1] on edges of triangles
303 SMESH_NodeXYZ myNode; //!< intersection node
304 bool myIsCollinear;//!< edges are collinear
306 IntPoint2D() : myIsCollinear( false ) {}
308 void InitLink( CutLink& link, int iFace, const std::vector< SMESH_NodeXYZ >& nodes ) const
310 link.Set( nodes[ myEdgeInd[ iFace ] ].Node(),
311 nodes[( myEdgeInd[ iFace ] + 1 ) % nodes.size() ].Node(),
313 link.myIntNode = myNode;
315 const SMDS_MeshNode* Node() const { return myNode.Node(); }
317 struct IntPoint2DCompare
320 IntPoint2DCompare( int iFace=0 ): myI( iFace ) {}
321 bool operator() ( const IntPoint2D* ip1, const IntPoint2D* ip2 ) const
323 return ip1->myU[ myI ] < ip2->myU[ myI ];
325 bool operator() ( const IntPoint2D& ip1, const IntPoint2D& ip2 ) const
327 return ip1.myU[ myI ] < ip2.myU[ myI ];
330 typedef boost::container::flat_set< IntPoint2D, IntPoint2DCompare > TIntPointSet;
331 typedef boost::container::flat_set< IntPoint2D*, IntPoint2DCompare > TIntPointPtrSet;
333 //--------------------------------------------------------------------------------
335 * \brief Face used to find translated position of the node
339 const SMDS_MeshElement* myFace;
340 SMESH_TNodeXYZ myNode1; //!< nodes neighboring another node of myFace
341 SMESH_TNodeXYZ myNode2;
342 const gp_XYZ* myNorm;
343 bool myNodeRightOrder;
344 void operator=(const SMDS_MeshElement* f) { myFace = f; }
345 const SMDS_MeshElement* operator->() { return myFace; }
346 void SetNodes( int i0, int i1 ) //!< set myNode's
348 myNode1.Set( myFace->GetNode( i1 ));
349 int i2 = ( i0 - 1 + myFace->NbCornerNodes() ) % myFace->NbCornerNodes();
351 i2 = ( i0 + 1 ) % myFace->NbCornerNodes();
352 myNode2.Set( myFace->GetNode( i2 ));
353 myNodeRightOrder = ( Abs( i2-i1 ) == 1 ) ? i2 > i1 : i2 < i1;
355 void SetOldNodes( const SMDS_Mesh& theSrcMesh )
357 myNode1.Set( theSrcMesh.FindNode( myNode1->GetID() ));
358 myNode2.Set( theSrcMesh.FindNode( myNode2->GetID() ));
360 bool SetNormal( const std::vector< gp_XYZ >& faceNormals )
362 myNorm = & faceNormals[ myFace->GetID() ];
363 return ( myNorm->SquareModulus() > gp::Resolution() * gp::Resolution() );
365 const gp_XYZ& Norm() const { return *myNorm; }
368 //--------------------------------------------------------------------------------
370 * \brief Offset plane used to find translated position of the node
377 gp_Lin myLines[2]; //!< line of intersection with neighbor OffsetPlane's
381 void Init( const gp_XYZ& node, Face& tria, double offset )
385 myPln = gp_Pln( node + tria.Norm() * offset, tria.Norm() );
386 myIsLineOk[0] = myIsLineOk[1] = false;
387 myWeight[0] = myWeight[1] = 0;
389 bool ComputeIntersectionLine( OffsetPlane& pln );
390 void SetSkewLine( const gp_Lin& line );
391 gp_XYZ GetCommonPoint( int & nbOkPoints, double& sumWeight );
392 gp_XYZ ProjectNodeOnLine( int & nbOkPoints );
393 double Weight() const { return myWeight[0] + myWeight[1]; }
396 //================================================================================
398 * \brief Set the second line
400 //================================================================================
402 void OffsetPlane::SetSkewLine( const gp_Lin& line )
405 gp_XYZ n = myLines[0].Direction().XYZ() ^ myLines[1].Direction().XYZ();
406 if (( myIsLineOk[1] = n.SquareModulus() > gp::Resolution() ))
407 myPln = gp_Pln( myPln.Location(), n );
410 //================================================================================
412 * \brief Project myNode on myLine[0]
414 //================================================================================
416 gp_XYZ OffsetPlane::ProjectNodeOnLine( int & nbOkPoints )
418 gp_XYZ p = gp::Origin().XYZ();
421 gp_Vec l2n( myLines[0].Location(), myNode );
422 double u = l2n * myLines[0].Direction();
423 p = myLines[0].Location().XYZ() + u * myLines[0].Direction().XYZ();
429 //================================================================================
431 * \brief Computes intersection point of myLines
433 //================================================================================
435 gp_XYZ OffsetPlane::GetCommonPoint( int & nbOkPoints, double& sumWeight )
437 if ( !myIsLineOk[0] || !myIsLineOk[1] )
439 // sumWeight += myWeight[0];
440 // return ProjectNodeOnLine( nbOkPoints ) * myWeight[0];
441 return gp::Origin().XYZ();
446 gp_Vec lPerp0 = myLines[0].Direction().XYZ() ^ myPln.Axis().Direction().XYZ();
447 double dot01 = lPerp0 * myLines[1].Direction().XYZ();
448 if ( Abs( dot01 ) > 0.05 )
450 gp_Vec l0l1 = myLines[1].Location().XYZ() - myLines[0].Location().XYZ();
451 double u1 = - ( lPerp0 * l0l1 ) / dot01;
452 p = ( myLines[1].Location().XYZ() + myLines[1].Direction().XYZ() * u1 );
456 gp_Vec lv0( myLines[0].Location(), myNode), lv1(myLines[1].Location(), myNode );
457 double dot0( lv0 * myLines[0].Direction() ), dot1( lv1 * myLines[1].Direction() );
458 p = 0.5 * ( myLines[0].Location().XYZ() + myLines[0].Direction().XYZ() * dot0 );
459 p += 0.5 * ( myLines[1].Location().XYZ() + myLines[1].Direction().XYZ() * dot1 );
462 sumWeight += Weight();
468 //================================================================================
470 * \brief Compute line of intersection of 2 planes
472 //================================================================================
474 bool OffsetPlane::ComputeIntersectionLine( OffsetPlane& theNextPln )
476 const gp_XYZ& n1 = myFace->Norm();
477 const gp_XYZ& n2 = theNextPln.myFace->Norm();
479 gp_XYZ lineDir = n1 ^ n2;
482 double x = Abs( lineDir.X() );
483 double y = Abs( lineDir.Y() );
484 double z = Abs( lineDir.Z() );
486 int cooMax; // max coordinate
488 if (x > z) cooMax = 1;
492 if (y > z) cooMax = 2;
497 if ( Abs( lineDir.Coord( cooMax )) < 0.05 )
499 // parallel planes - intersection is an offset of the common edge
500 linePos = 0.5 * ( myPln.Location().XYZ() + theNextPln.myPln.Location().XYZ() );
501 lineDir = myNode - myFace->myNode2;
507 // the constants in the 2 plane equations
508 double d1 = - ( n1 * myPln.Location().XYZ() );
509 double d2 = - ( n2 * theNextPln.myPln.Location().XYZ() );
514 linePos.SetY(( d2*n1.Z() - d1*n2.Z()) / lineDir.X() );
515 linePos.SetZ(( d1*n2.Y() - d2*n1.Y()) / lineDir.X() );
518 linePos.SetX(( d1*n2.Z() - d2*n1.Z()) / lineDir.Y() );
520 linePos.SetZ(( d2*n1.X() - d1*n2.X()) / lineDir.Y() );
523 linePos.SetX(( d2*n1.Y() - d1*n2.Y()) / lineDir.Z() );
524 linePos.SetY(( d1*n2.X() - d2*n1.X()) / lineDir.Z() );
527 myWeight[0] = lineDir.SquareModulus();
529 myWeight[0] = 2. - myWeight[0];
531 myLines [ 0 ].SetDirection( lineDir );
532 myLines [ 0 ].SetLocation ( linePos );
533 myIsLineOk[ 0 ] = ok;
535 theNextPln.myLines [ 1 ] = myLines[ 0 ];
536 theNextPln.myIsLineOk[ 1 ] = ok;
537 theNextPln.myWeight [ 1 ] = myWeight[ 0 ];
542 //================================================================================
544 * \brief Return a translated position of a node
545 * \param [in] new2OldNodes - new and old nodes
546 * \param [in] faceNormals - normals to input faces
547 * \param [in] theSrcMesh - initial mesh
548 * \param [in] theNewPos - a computed normal
549 * \return bool - true if theNewPos is computed
551 //================================================================================
553 bool getTranslatedPosition( const SMDS_MeshNode* theNewNode,
554 const double theOffset,
556 const double theSign,
557 const std::vector< gp_XYZ >& theFaceNormals,
558 SMDS_Mesh& theSrcMesh,
561 bool useOneNormal = true;
563 // check if theNewNode needs an average position, i.e. theNewNode is convex
564 // SMDS_ElemIteratorPtr faceIt = theNewNode->GetInverseElementIterator();
565 // const SMDS_MeshElement* f0 = faceIt->next();
566 // const gp_XYZ& norm0 = theFaceNormals[ f0->GetID() ];
567 // const SMESH_NodeXYZ nodePos = theNewNode;
568 // while ( faceIt->more() )
570 // const SMDS_MeshElement* f = faceIt->next();
571 // const int nodeInd = f->GetNodeIndex( theNewNode );
572 // SMESH_NodeXYZ nodePos2 = f->GetWrapNode( nodeInd + 1 );
574 // const gp_XYZ nnDir = ( nodePos2 - nodePos ).Normalized();
579 // const double dot = norm0 * nnDir;
584 // get faces surrounding theNewNode and sort them
585 Face faces[ theMaxNbFaces ];
586 SMDS_ElemIteratorPtr faceIt = theNewNode->GetInverseElementIterator();
587 faces[0] = faceIt->next();
588 while ( !faces[0].SetNormal( theFaceNormals ) && faceIt->more() )
589 faces[0] = faceIt->next();
590 int i0 = faces[0]->GetNodeIndex( theNewNode );
591 int i1 = ( i0 + 1 ) % faces[0]->NbCornerNodes();
592 faces[0].SetNodes( i0, i1 );
593 TIDSortedElemSet elemSet, avoidSet;
595 const SMDS_MeshElement* f;
596 for ( ; faceIt->more() && iFace < theMaxNbFaces; faceIt->next() )
598 avoidSet.insert( faces[ iFace ].myFace );
599 f = SMESH_MeshAlgos::FindFaceInSet( theNewNode, faces[ iFace ].myNode2.Node(),
600 elemSet, avoidSet, &i0, &i1 );
603 std::reverse( &faces[0], &faces[0] + iFace + 1 );
604 for ( int i = 0; i <= iFace; ++i )
606 std::swap( faces[i].myNode1, faces[i].myNode2 );
607 faces[i].myNodeRightOrder = !faces[i].myNodeRightOrder;
609 f = SMESH_MeshAlgos::FindFaceInSet( theNewNode, faces[ iFace ].myNode2.Node(),
610 elemSet, avoidSet, &i0, &i1 );
614 faces[ ++iFace ] = f;
615 faces[ iFace ].SetNodes( i0, i1 );
616 faces[ iFace ].SetNormal( theFaceNormals );
618 int nbFaces = iFace + 1;
620 theNewPos.SetCoord( 0, 0, 0 );
621 gp_XYZ oldXYZ = SMESH_NodeXYZ( theNewNode );
623 // check if all faces are co-planar
624 bool isPlanar = true;
625 const double tol = 1e-2;
626 for ( int i = 1; i < nbFaces && isPlanar; ++i )
627 isPlanar = ( faces[i].Norm() - faces[i-1].Norm() ).SquareModulus() < tol*tol;
631 theNewPos = oldXYZ + faces[0].Norm() * theOffset;
635 // prepare OffsetPlane's
636 OffsetPlane pln[ theMaxNbFaces ];
637 for ( int i = 0; i < nbFaces; ++i )
639 faces[i].SetOldNodes( theSrcMesh );
640 pln[i].Init( oldXYZ, faces[i], theOffset );
642 // intersect neighboring OffsetPlane's
644 for ( int i = 1; i < nbFaces; ++i )
645 nbOkPoints += pln[ i-1 ].ComputeIntersectionLine( pln[ i ]);
646 nbOkPoints += pln[ nbFaces-1 ].ComputeIntersectionLine( pln[ 0 ]);
648 // move intersection lines to over parallel planes
649 if ( nbOkPoints > 1 )
650 for ( int i = 0; i < nbFaces; ++i )
651 if ( pln[i].myIsLineOk[0] && !pln[i].myIsLineOk[1] )
652 for ( int j = 1; j < nbFaces && !pln[i].myIsLineOk[1]; ++j )
654 int i2 = ( i + j ) % nbFaces;
655 if ( pln[i2].myIsLineOk[0] )
656 pln[i].SetSkewLine( pln[i2].myLines[0] );
659 // get the translated position
661 double sumWegith = 0;
662 const double minWeight = Sin( 30 * M_PI / 180. ) * Sin( 30 * M_PI / 180. );
663 for ( int i = 0; i < nbFaces; ++i )
664 if ( pln[ i ].Weight() > minWeight )
665 theNewPos += pln[ i ].GetCommonPoint( nbOkPoints, sumWegith );
667 if ( nbOkPoints == 0 )
669 // there is only one feature edge;
670 // find the theNewPos by projecting oldXYZ to any intersection line
671 for ( int i = 0; i < nbFaces; ++i )
672 theNewPos += pln[ i ].ProjectNodeOnLine( nbOkPoints );
674 if ( nbOkPoints == 0 )
676 theNewPos = oldXYZ + faces[0].Norm() * theOffset;
679 sumWegith = nbOkPoints;
681 theNewPos /= sumWegith;
684 // mark theNewNode if it is concave
685 useOneNormal = false;
686 gp_Vec moveVec( oldXYZ, theNewPos );
687 for ( int i = 0, iPrev = nbFaces-1; i < nbFaces; iPrev = i++ )
689 gp_Vec nodeVec( oldXYZ, faces[ i ].myNode1 );
690 double u = ( moveVec * nodeVec ) / nodeVec.SquareMagnitude();
691 if ( u > 0.5 ) // param [0,1] on nodeVec
693 theNewNode->setIsMarked( true );
697 gp_XYZ inFaceVec = faces[ i ].Norm() ^ nodeVec.XYZ();
698 double dot = inFaceVec * faces[ iPrev ].Norm();
699 if ( !faces[ i ].myNodeRightOrder )
701 if ( dot * theSign < 0 )
703 gp_XYZ p1 = oldXYZ + faces[ i ].Norm() * theOffset;
704 gp_XYZ p2 = oldXYZ + faces[ iPrev ].Norm() * theOffset;
705 useOneNormal = ( p1 - p2 ).SquareModulus() > 1e-12;
708 if ( useOneNormal && theNewNode->isMarked() )
715 //================================================================================
717 * \brief Remove small faces
719 //================================================================================
721 void removeSmallFaces( SMDS_Mesh* theMesh,
722 SMESH_MeshAlgos::TElemIntPairVec& theNew2OldFaces,
723 const double theTol2 )
725 std::vector< SMESH_NodeXYZ > points(3);
726 std::vector< const SMDS_MeshNode* > nodes(3);
727 for ( SMDS_ElemIteratorPtr faceIt = theMesh->elementsIterator(); faceIt->more(); )
729 const SMDS_MeshElement* face = faceIt->next();
730 points.assign( face->begin_nodes(), face->end_nodes() );
732 SMESH_NodeXYZ* prevN = & points.back();
733 for ( size_t i = 0; i < points.size(); ++i )
735 double dist2 = ( *prevN - points[ i ]).SquareModulus();
736 if ( dist2 < theTol2 )
738 const SMDS_MeshNode* nToRemove =
739 (*prevN)->GetID() > points[ i ]->GetID() ? prevN->Node() : points[ i ].Node();
740 const SMDS_MeshNode* nToKeep =
741 nToRemove == points[ i ].Node() ? prevN->Node() : points[ i ].Node();
742 for ( SMDS_ElemIteratorPtr fIt = nToRemove->GetInverseElementIterator(); fIt->more(); )
744 const SMDS_MeshElement* f = fIt->next();
747 nodes.assign( f->begin_nodes(), f->end_nodes() );
748 nodes[ f->GetNodeIndex( nToRemove )] = nToKeep;
749 theMesh->ChangeElementNodes( f, &nodes[0], nodes.size() );
751 theNew2OldFaces[ face->GetID() ].first = 0;
752 theMesh->RemoveFreeElement( face );
755 prevN = & points[ i ];
764 namespace SMESH_MeshAlgos
766 //--------------------------------------------------------------------------------
768 * \brief Intersect faces of a mesh
770 struct Intersector::Algo
774 const std::vector< gp_XYZ >& myNormals;
775 TCutLinkMap myCutLinks; //!< assure sharing of new nodes
776 TCutFaceMap myCutFaces;
777 TNNMap myRemove2KeepNodes; //!< node merge map
779 // data to intersect 2 faces
780 const SMDS_MeshElement* myFace1;
781 const SMDS_MeshElement* myFace2;
782 std::vector< SMESH_NodeXYZ > myNodes1, myNodes2;
783 std::vector< double > myDist1, myDist2;
784 int myInd1, myInd2; // coordinate indices on an axis-aligned plane
785 int myNbOnPlane1, myNbOnPlane2;
786 TIntPointSet myIntPointSet;
788 Algo( SMDS_Mesh* mesh, double tol, const std::vector< gp_XYZ >& normals )
792 //myEps( Sqrt( std::numeric_limits<double>::min() )),
793 //myEps( gp::Resolution() ),
796 void Cut( const SMDS_MeshElement* face1,
797 const SMDS_MeshElement* face2,
798 const int nbCommonNodes );
799 void Cut( const SMDS_MeshElement* face,
800 SMESH_NodeXYZ& lineEnd1,
802 SMESH_NodeXYZ& lineEnd2,
804 void MakeNewFaces( TElemIntPairVec& theNew2OldFaces,
805 TNodeIntPairVec& theNew2OldNodes,
806 const double theSign,
807 const bool theOptimize );
809 void IntersectNewEdges( const CutFace& theCFace );
813 bool isPlaneIntersected( const gp_XYZ& n2,
815 const std::vector< SMESH_NodeXYZ >& nodes1,
816 std::vector< double > & dist1,
819 void computeIntervals( const std::vector< SMESH_NodeXYZ >& nodes,
820 const std::vector< double >& dist,
826 void addLink ( CutLink& link );
827 bool findLink( CutLink& link );
828 bool coincide( const gp_XYZ& p1, const gp_XYZ& p2, const double tol ) const
830 return ( p1 - p2 ).SquareModulus() < tol * tol;
832 gp_XY p2D( const gp_XYZ& p ) const { return gp_XY( p.Coord( myInd1 ), p.Coord( myInd2 )); }
834 void intersectLink( const std::vector< SMESH_NodeXYZ >& nodes1,
835 const std::vector< double > & dist1,
837 const SMDS_MeshElement* face2,
839 void findIntPointOnPlane( const std::vector< SMESH_NodeXYZ >& nodes,
840 const std::vector< double > & dist,
842 void replaceIntNode( const SMDS_MeshNode* nToKeep, const SMDS_MeshNode* nToRemove );
843 void computeIntPoint( const double u1,
848 const SMDS_MeshNode* & node1,
849 const SMDS_MeshNode* & node2);
850 void cutCollinearLink( const int iNotOnPlane1,
851 const std::vector< SMESH_NodeXYZ >& nodes1,
852 const SMDS_MeshElement* face2,
853 const CutLink& link1,
854 const CutLink& link2);
855 void setPlaneIndices( const gp_XYZ& planeNorm );
856 bool intersectEdgeEdge( const gp_XY s1p0, const gp_XY s1p1,
857 const gp_XY s2p0, const gp_XY s2p1,
858 double & t1, double & t2,
859 bool & isCollinear );
860 bool intersectEdgeEdge( int iE1, int iE2, IntPoint2D& intPoint );
861 bool isPointInTriangle( const gp_XYZ& p, const std::vector< SMESH_NodeXYZ >& nodes );
862 const SMDS_MeshNode* createNode( const gp_XYZ& p );
865 //================================================================================
867 * \brief Return coordinate index with maximal abs value
869 //================================================================================
871 int MaxIndex( const gp_XYZ& x )
873 int iMaxCoo = ( Abs( x.X()) < Abs( x.Y() )) + 1;
874 if ( Abs( x.Coord( iMaxCoo )) < Abs( x.Z() ))
878 //================================================================================
880 * \brief Store a CutLink
882 //================================================================================
884 const SMDS_MeshNode* Intersector::Algo::createNode( const gp_XYZ& p )
886 const SMDS_MeshNode* n = myMesh->AddNode( p.X(), p.Y(), p.Z() );
887 n->setIsMarked( true ); // cut nodes are marked
891 //================================================================================
893 * \brief Store a CutLink
895 //================================================================================
897 void Intersector::Algo::addLink( CutLink& link )
900 const CutLink* added = & myCutLinks.Added( link );
901 while ( added->myIntNode.Node() != link.myIntNode.Node() )
903 if ( !added->myIntNode )
905 added->myIntNode = link.myIntNode;
911 added = & myCutLinks.Added( link );
917 //================================================================================
919 * \brief Find a CutLink with an intersection point coincident with that of a given link
921 //================================================================================
923 bool Intersector::Algo::findLink( CutLink& link )
926 while ( myCutLinks.Contains( link ))
928 const CutLink* added = & myCutLinks.Added( link );
929 if ( !!added->myIntNode && coincide( added->myIntNode, link.myIntNode, myTol ))
931 link.myIntNode = added->myIntNode;
939 //================================================================================
941 * \brief Check if a triangle intersects the plane of another triangle
942 * \param [in] nodes1 - nodes of triangle 1
943 * \param [in] n2 - normal of triangle 2
944 * \param [in] d2 - a constant of the plane equation 2
945 * \param [out] dist1 - distance of nodes1 from the plane 2
946 * \param [out] nbOnPlane - number of nodes1 lying on the plane 2
947 * \return bool - true if the triangle intersects the plane 2
949 //================================================================================
951 bool Intersector::Algo::isPlaneIntersected( const gp_XYZ& n2,
953 const std::vector< SMESH_NodeXYZ >& nodes1,
954 std::vector< double > & dist1,
958 iNotOnPlane1 = nbOnPlane1 = 0;
959 dist1.resize( nodes1.size() );
960 for ( size_t i = 0; i < nodes1.size(); ++i )
962 dist1[i] = n2 * nodes1[i] + d2;
963 if ( Abs( dist1[i] ) < myTol )
973 if ( nbOnPlane1 == 0 )
974 for ( size_t i = 0; i < nodes1.size(); ++i )
975 if ( dist1[iNotOnPlane1] * dist1[i] < 0 )
981 //================================================================================
983 * \brief Compute parameters on the plane intersection line of intersections
984 * of edges of a triangle
985 * \param [in] nodes - triangle nodes
986 * \param [in] dist - distance of triangle nodes from the plane of another triangle
987 * \param [in] nbOnPln - number of nodes lying on the plane of another triangle
988 * \param [in] iMaxCoo - index of coordinate of max component of the plane intersection line
989 * \param [out] u - two computed parameters on the plane intersection line
990 * \param [out] iE - indices of intersected edges
992 //================================================================================
994 void Intersector::Algo::computeIntervals( const std::vector< SMESH_NodeXYZ >& nodes,
995 const std::vector< double >& dist,
1003 u[0] = u[1] = 1e+100;
1008 if ( nbOnPln == 1 && ( dist[i1] == 0. || dist[i2] == 0 ))
1010 int i = dist[i1] == 0 ? i1 : i2;
1011 u [ 1 ] = nodes[ i ].Coord( iMaxCoo );
1015 for ( ; i2 < 3 && nb < 2; i1 = i2++ )
1017 double dd = dist[i1] - dist[i2];
1018 if ( dd != 0. && dist[i2] * dist[i1] <= 0. )
1020 double x1 = nodes[i1].Coord( iMaxCoo );
1021 double x2 = nodes[i2].Coord( iMaxCoo );
1022 u [ nb ] = x1 + ( x2 - x1 ) * dist[i1] / dd;
1029 std::swap( u [0], u [1] );
1030 std::swap( iE[0], iE[1] );
1034 //================================================================================
1036 * \brief Try to find an intersection node on a link collinear with the plane intersection line
1038 //================================================================================
1040 void Intersector::Algo::findIntPointOnPlane( const std::vector< SMESH_NodeXYZ >& nodes,
1041 const std::vector< double > & dist,
1044 int i1 = ( dist[0] == 0 ? 0 : 1 ), i2 = ( dist[2] == 0 ? 2 : 1 );
1045 CutLink link2 = link;
1046 link2.Set( nodes[i1].Node(), nodes[i2].Node(), 0 );
1047 if ( findLink( link2 ))
1048 link.myIntNode = link2.myIntNode;
1051 //================================================================================
1053 * \brief Compute intersection point of a link1 with a face2
1055 //================================================================================
1057 void Intersector::Algo::intersectLink( const std::vector< SMESH_NodeXYZ >& nodes1,
1058 const std::vector< double > & dist1,
1060 const SMDS_MeshElement* face2,
1063 const int iEdge2 = ( iEdge1 + 1 ) % nodes1.size();
1064 const SMESH_NodeXYZ& p1 = nodes1[ iEdge1 ];
1065 const SMESH_NodeXYZ& p2 = nodes1[ iEdge2 ];
1067 link1.Set( p1.Node(), p2.Node(), face2 );
1068 const CutLink* link = & myCutLinks.Added( link1 );
1069 if ( !link->IntNode() )
1071 if ( dist1[ iEdge1 ] == 0. ) link1.myIntNode = p1;
1072 else if ( dist1[ iEdge2 ] == 0. ) link1.myIntNode = p2;
1075 gp_XYZ p = p1 + ( p2 - p1 ) * dist1[ iEdge1 ] / ( dist1[ iEdge1 ] - dist1[ iEdge2 ]);
1076 (gp_XYZ&)link1.myIntNode = p;
1081 gp_XYZ p = p1 + ( p2 - p1 ) * dist1[ iEdge1 ] / ( dist1[ iEdge1 ] - dist1[ iEdge2 ]);
1082 while ( link->IntNode() )
1084 if ( coincide( p, link->myIntNode, myTol ))
1086 link1.myIntNode = link->myIntNode;
1090 link = & myCutLinks.Added( link1 );
1092 if ( !link1.IntNode() )
1094 if ( dist1[ iEdge1 ] == 0. ) link1.myIntNode = p1;
1095 else if ( dist1[ iEdge2 ] == 0. ) link1.myIntNode = p2;
1096 else (gp_XYZ&)link1.myIntNode = p;
1101 //================================================================================
1103 * \brief Store node replacement in myCutFaces
1105 //================================================================================
1107 void Intersector::Algo::replaceIntNode( const SMDS_MeshNode* nToKeep,
1108 const SMDS_MeshNode* nToRemove )
1110 if ( nToKeep == nToRemove )
1112 if ( nToRemove->GetID() < nToKeep->GetID() ) // keep node with lower ID
1113 myRemove2KeepNodes.Bind( nToKeep, nToRemove );
1115 myRemove2KeepNodes.Bind( nToRemove, nToKeep );
1118 //================================================================================
1120 * \brief Compute intersection point on a link of either of faces by choosing
1121 * a link whose parameter on the intersection line in maximal
1122 * \param [in] u1 - parameter on the intersection line of link iE1 of myFace1
1123 * \param [in] u2 - parameter on the intersection line of link iE2 of myFace2
1124 * \param [in] iE1 - index of a link myFace1
1125 * \param [in] iE2 - index of a link myFace2
1126 * \param [out] link - CutLink storing the intersection point
1127 * \param [out] node1 - a node of the 2nd link if two links intersect
1128 * \param [out] node2 - a node of the 2nd link if two links intersect
1130 //================================================================================
1132 void Intersector::Algo::computeIntPoint( const double u1,
1137 const SMDS_MeshNode* & node1,
1138 const SMDS_MeshNode* & node2)
1140 if ( u1 > u2 + myTol )
1142 intersectLink( myNodes1, myDist1, iE1, myFace2, link );
1144 if ( myNbOnPlane2 == 2 )
1145 findIntPointOnPlane( myNodes2, myDist2, link );
1147 else if ( u2 > u1 + myTol )
1149 intersectLink( myNodes2, myDist2, iE2, myFace1, link );
1151 if ( myNbOnPlane1 == 2 )
1152 findIntPointOnPlane( myNodes1, myDist1, link );
1154 else // edges of two faces intersect the line at the same point
1157 intersectLink( myNodes1, myDist1, iE1, myFace2, link );
1158 intersectLink( myNodes2, myDist2, iE2, myFace1, link2 );
1159 node1 = link2.Node1();
1160 node2 = link2.Node2();
1162 if ( !link.IntNode() && link2.IntNode() )
1163 link.myIntNode = link2.myIntNode;
1165 else if ( !link.IntNode() && !link2.IntNode() )
1166 (gp_XYZ&)link.myIntNode = 0.5 * ( link.myIntNode + link2.myIntNode );
1168 else if ( link.IntNode() && link2.IntNode() )
1169 replaceIntNode( link.IntNode(), link2.IntNode() );
1173 //================================================================================
1175 * \brief Add intersections to a link collinear with the intersection line
1177 //================================================================================
1179 void Intersector::Algo::cutCollinearLink( const int iNotOnPlane1,
1180 const std::vector< SMESH_NodeXYZ >& nodes1,
1181 const SMDS_MeshElement* face2,
1182 const CutLink& link1,
1183 const CutLink& link2)
1186 int iN1 = ( iNotOnPlane1 + 1 ) % 3;
1187 int iN2 = ( iNotOnPlane1 + 2 ) % 3;
1188 CutLink link( nodes1[ iN1 ].Node(), nodes1[ iN2 ].Node(), face2 );
1189 if ( link1.myFace != face2 )
1191 link.myIntNode = link1.myIntNode;
1194 if ( link2.myFace != face2 )
1196 link.myIntNode = link2.myIntNode;
1201 //================================================================================
1203 * \brief Choose indices on an axis-aligned plane
1205 //================================================================================
1207 void Intersector::Algo::setPlaneIndices( const gp_XYZ& planeNorm )
1209 switch ( MaxIndex( planeNorm )) {
1210 case 1: myInd1 = 2; myInd2 = 3; break;
1211 case 2: myInd1 = 3; myInd2 = 1; break;
1212 case 3: myInd1 = 1; myInd2 = 2; break;
1216 //================================================================================
1218 * \brief Intersect two faces
1220 //================================================================================
1222 void Intersector::Algo::Cut( const SMDS_MeshElement* face1,
1223 const SMDS_MeshElement* face2,
1224 const int nbCommonNodes)
1228 myNodes1.assign( face1->begin_nodes(), face1->end_nodes() );
1229 myNodes2.assign( face2->begin_nodes(), face2->end_nodes() );
1231 const gp_XYZ& n1 = myNormals[ face1->GetID() ];
1232 const gp_XYZ& n2 = myNormals[ face2->GetID() ];
1234 // check if triangles intersect
1235 int iNotOnPlane1, iNotOnPlane2;
1236 const double d2 = -( n2 * myNodes2[0]);
1237 if ( !isPlaneIntersected( n2, d2, myNodes1, myDist1, myNbOnPlane1, iNotOnPlane1 ))
1239 const double d1 = -( n1 * myNodes1[0]);
1240 if ( !isPlaneIntersected( n1, d1, myNodes2, myDist2, myNbOnPlane2, iNotOnPlane2 ))
1243 if ( myNbOnPlane1 == 3 || myNbOnPlane2 == 3 )// triangles are co-planar
1245 setPlaneIndices( myNbOnPlane1 == 3 ? n2 : n1 ); // choose indices on an axis-aligned plane
1248 else if ( nbCommonNodes < 2 ) // triangle planes intersect
1250 gp_XYZ lineDir = n1 ^ n2; // intersection line
1252 // check if intervals of intersections of triangles with lineDir overlap
1254 double u1[2], u2 [2]; // parameters on lineDir of edge intersection points { minU, maxU }
1255 int iE1[2], iE2[2]; // indices of edges
1256 int iMaxCoo = MaxIndex( lineDir );
1257 computeIntervals( myNodes1, myDist1, myNbOnPlane1, iMaxCoo, u1, iE1 );
1258 computeIntervals( myNodes2, myDist2, myNbOnPlane2, iMaxCoo, u2, iE2 );
1259 if ( u1[1] < u2[0] - myTol || u2[1] < u1[0] - myTol )
1260 return; // intervals do not overlap
1262 // make intersection nodes
1264 const SMDS_MeshNode *l1n1, *l1n2, *l2n1, *l2n2;
1265 CutLink link1; // intersection with smaller u on lineDir
1266 computeIntPoint( u1[0], u2[0], iE1[0], iE2[0], link1, l1n1, l1n2 );
1267 CutLink link2; // intersection with larger u on lineDir
1268 computeIntPoint( -u1[1], -u2[1], iE1[1], iE2[1], link2, l2n1, l2n2 );
1270 const CutFace& cf1 = myCutFaces.Added( CutFace( face1 ));
1271 const CutFace& cf2 = myCutFaces.Added( CutFace( face2 ));
1273 if ( coincide( link1.myIntNode, link2.myIntNode, myTol ))
1275 // intersection is a point
1276 if ( link1.IntNode() && link2.IntNode() )
1277 replaceIntNode( link1.IntNode(), link2.IntNode() );
1279 CutLink* link = link2.IntNode() ? &link2 : &link1;
1280 if ( !link->IntNode() )
1282 gp_XYZ p = 0.5 * ( link1.myIntNode + link2.myIntNode );
1283 link->myIntNode.Set( createNode( p ));
1285 if ( !link1.IntNode() ) link1.myIntNode = link2.myIntNode;
1286 if ( !link2.IntNode() ) link2.myIntNode = link1.myIntNode;
1288 cf1.AddPoint( link1, link2, myTol );
1289 if ( l1n1 ) link1.Set( l1n1, l1n2, face2 );
1290 if ( l2n1 ) link2.Set( l2n1, l2n2, face2 );
1291 cf2.AddPoint( link1, link2, myTol );
1295 // intersection is a line segment
1296 if ( !link1.IntNode() )
1297 link1.myIntNode.Set( createNode( link1.myIntNode ));
1298 if ( !link2.IntNode() )
1299 link2.myIntNode.Set( createNode( link2.myIntNode ));
1301 cf1.AddEdge( link1, link2, face2, myNbOnPlane1, iNotOnPlane1 );
1302 if ( l1n1 ) link1.Set( l1n1, l1n2, face2 );
1303 if ( l2n1 ) link2.Set( l2n1, l2n2, face2 );
1304 cf2.AddEdge( link1, link2, face1, myNbOnPlane2, iNotOnPlane2 );
1306 // add intersections to a link collinear with the intersection line
1307 if ( myNbOnPlane1 == 2 && ( link1.myFace != face2 || link2.myFace != face2 ))
1308 cutCollinearLink( iNotOnPlane1, myNodes1, face2, link1, link2 );
1310 if ( myNbOnPlane2 == 2 && ( link1.myFace != face1 || link2.myFace != face1 ))
1311 cutCollinearLink( iNotOnPlane2, myNodes2, face1, link1, link2 );
1317 } // non co-planar case
1322 //================================================================================
1324 * \brief Store a face cut by a line given by its ends
1325 * accompanied by indices of intersected face edges.
1326 * Edge index is <0 if a line end is inside the face.
1327 * \param [in] face - a face to cut
1328 * \param [inout] lineEnd1 - line end coordinates + optional node existing at this point
1329 * \param [in] edgeIndex1 - index of face edge cut by lineEnd1
1330 * \param [inout] lineEnd2 - line end coordinates + optional node existing at this point
1331 * \param [in] edgeIndex2 - index of face edge cut by lineEnd2
1333 //================================================================================
1335 void Intersector::Algo::Cut( const SMDS_MeshElement* face,
1336 SMESH_NodeXYZ& lineEnd1,
1338 SMESH_NodeXYZ& lineEnd2,
1341 if ( lineEnd1.Node() && lineEnd2.Node() &&
1342 face->GetNodeIndex( lineEnd1.Node() ) >= 0 &&
1343 face->GetNodeIndex( lineEnd2.Node() ) >= 0 )
1344 return; // intersection at a face node or edge
1346 if ((int) myNormals.size() <= face->GetID() )
1347 const_cast< std::vector< gp_XYZ >& >( myNormals ).resize( face->GetID() + 1 );
1349 const CutFace& cf = myCutFaces.Added( CutFace( face ));
1352 // look for intersection nodes coincident with line ends
1354 for ( int is2nd = 0; is2nd < 2; ++is2nd )
1356 SMESH_NodeXYZ& lineEnd = is2nd ? lineEnd2 : lineEnd1;
1357 int edgeIndex = is2nd ? edgeIndex2 : edgeIndex1;
1358 CutLink & link = links[ is2nd ];
1360 link.myIntNode = lineEnd;
1362 for ( size_t i = ( edgeIndex < 0 ? 3 : 0 ); i < cf.myLinks.size(); ++i )
1363 if ( coincide( lineEnd, SMESH_NodeXYZ( cf.myLinks[i].myNode1 ), myTol ))
1365 link.myIntNode = cf.myLinks[i].myNode1;
1369 if ( edgeIndex >= 0 )
1371 link.Set( face->GetNode ( edgeIndex ),
1372 face->GetNodeWrap( edgeIndex + 1 ),
1377 if ( !link.myIntNode )
1378 link.myIntNode.Set( createNode( lineEnd ));
1380 lineEnd._node = link.IntNode();
1382 if ( link.myNode[0] )
1386 cf.AddEdge( links[0], links[1], /*face=*/0, /*nbOnPlane=*/0, /*iNotOnPlane=*/-1 );
1389 //================================================================================
1391 * \brief Intersect two 2D line segments
1393 //================================================================================
1395 bool Intersector::Algo::intersectEdgeEdge( const gp_XY s1p0, const gp_XY s1p1,
1396 const gp_XY s2p0, const gp_XY s2p1,
1397 double & t1, double & t2,
1398 bool & isCollinear )
1400 gp_XY u = s1p1 - s1p0;
1401 gp_XY v = s2p1 - s2p0;
1402 gp_XY w = s1p0 - s2p0;
1403 double perpDotUV = u * gp_XY( -v.Y(), v.X() );
1404 double perpDotVW = v * gp_XY( -w.Y(), w.X() );
1405 double perpDotUW = u * gp_XY( -w.Y(), w.X() );
1406 double u2 = u.SquareModulus();
1407 double v2 = v.SquareModulus();
1408 if ( u2 < myEps * myEps || v2 < myEps * myEps )
1410 if ( perpDotUV * perpDotUV / u2 / v2 < 1e-6 ) // cos ^ 2
1413 return false; // no need in collinear solution
1414 if ( perpDotUW * perpDotUW / u2 > myTol * myTol )
1415 return false; // parallel
1418 gp_XY w2 = s1p1 - s2p0;
1419 if ( Abs( v.X()) + Abs( u.X()) > Abs( v.Y()) + Abs( u.Y())) {
1420 t1 = w.X() / v.X(); // params on segment 2
1421 t2 = w2.X() / v.X();
1425 t2 = w2.Y() / v.Y();
1427 if ( Max( t1,t2 ) <= 0 || Min( t1,t2 ) >= 1 )
1428 return false; // no overlap
1431 isCollinear = false;
1433 t1 = perpDotVW / perpDotUV; // param on segment 1
1434 if ( t1 < 0. || t1 > 1. )
1435 return false; // intersection not within the segment
1437 t2 = perpDotUW / perpDotUV; // param on segment 2
1438 if ( t2 < 0. || t2 > 1. )
1439 return false; // intersection not within the segment
1444 //================================================================================
1446 * \brief Intersect two edges of co-planar triangles
1447 * \param [inout] iE1 - edge index of triangle 1
1448 * \param [inout] iE2 - edge index of triangle 2
1449 * \param [inout] intPoints - intersection points
1450 * \param [inout] nbIntPoints - nb of found intersection points
1452 //================================================================================
1454 bool Intersector::Algo::intersectEdgeEdge( int iE1, int iE2, IntPoint2D& intPoint )
1456 int i01 = iE1, i11 = ( iE1 + 1 ) % 3;
1457 int i02 = iE2, i12 = ( iE2 + 1 ) % 3;
1458 if (( !intPoint.myIsCollinear ) &&
1459 ( myNodes1[ i01 ] == myNodes2[ i02 ] ||
1460 myNodes1[ i01 ] == myNodes2[ i12 ] ||
1461 myNodes1[ i11 ] == myNodes2[ i02 ] ||
1462 myNodes1[ i11 ] == myNodes2[ i12 ] ))
1466 gp_XY s1p0 = p2D( myNodes1[ i01 ]);
1467 gp_XY s1p1 = p2D( myNodes1[ i11 ]);
1470 gp_XY s2p0 = p2D( myNodes2[ i02 ]);
1471 gp_XY s2p1 = p2D( myNodes2[ i12 ]);
1474 if ( !intersectEdgeEdge( s1p0,s1p1, s2p0,s2p1, t1, t2, intPoint.myIsCollinear ))
1477 intPoint.myEdgeInd[0] = iE1;
1478 intPoint.myEdgeInd[1] = iE2;
1479 intPoint.myU[0] = t1;
1480 intPoint.myU[1] = t2;
1481 (gp_XYZ&)intPoint.myNode = myNodes1[i01] * ( 1 - t1 ) + myNodes1[i11] * t1;
1483 if ( intPoint.myIsCollinear )
1486 // try to find existing node at intPoint.myNode
1488 if ( myNodes1[ i01 ] == myNodes2[ i02 ] ||
1489 myNodes1[ i01 ] == myNodes2[ i12 ] ||
1490 myNodes1[ i11 ] == myNodes2[ i02 ] ||
1491 myNodes1[ i11 ] == myNodes2[ i12 ] )
1494 const double coincTol = myTol * 1e-3;
1496 CutLink link1( myNodes1[i01].Node(), myNodes1[i11].Node(), myFace2 );
1497 CutLink link2( myNodes2[i02].Node(), myNodes2[i12].Node(), myFace1 );
1499 SMESH_NodeXYZ& n1 = myNodes1[ t1 < 0.5 ? i01 : i11 ];
1500 bool same1 = coincide( n1, intPoint.myNode, coincTol );
1503 link2.myIntNode = intPoint.myNode = n1;
1506 SMESH_NodeXYZ& n2 = myNodes2[ t2 < 0.5 ? i02 : i12 ];
1507 bool same2 = coincide( n2, intPoint.myNode, coincTol );
1510 link1.myIntNode = intPoint.myNode = n2;
1514 replaceIntNode( n1.Node(), n2.Node() );
1522 link1.myIntNode = intPoint.myNode;
1523 if ( findLink( link1 ))
1525 intPoint.myNode = link2.myIntNode = link1.myIntNode;
1530 link2.myIntNode = intPoint.myNode;
1531 if ( findLink( link2 ))
1533 intPoint.myNode = link1.myIntNode = link2.myIntNode;
1538 for ( int is2nd = 0; is2nd < 2; ++is2nd )
1540 const SMDS_MeshElement* f = is2nd ? myFace1 : myFace2;
1542 const CutFace& cf = myCutFaces.Added( CutFace( is2nd ? myFace2 : myFace1 ));
1543 for ( size_t i = 0; i < cf.myLinks.size(); ++i )
1544 if ( cf.myLinks[i].myFace == f &&
1545 //cf.myLinks[i].myIndex != EdgePart::_COPLANAR &&
1546 coincide( intPoint.myNode, SMESH_NodeXYZ( cf.myLinks[i].myNode1 ), coincTol ))
1548 intPoint.myNode.Set( cf.myLinks[i].myNode1 );
1555 intPoint.myNode._node = createNode( intPoint.myNode );
1556 link1.myIntNode = link2.myIntNode = intPoint.myNode;
1564 //================================================================================
1566 * \brief Check if a point is contained in a triangle
1568 //================================================================================
1570 bool Intersector::Algo::isPointInTriangle( const gp_XYZ& p, const std::vector< SMESH_NodeXYZ >& nodes )
1573 SMESH_MeshAlgos::GetBarycentricCoords( p2D( p ),
1574 p2D( nodes[0] ), p2D( nodes[1] ), p2D( nodes[2] ),
1576 //return ( 0. < bc1 && 0. < bc2 && bc1 + bc2 < 1. );
1577 return ( myTol < bc1 && myTol < bc2 && bc1 + bc2 + myTol < 1. );
1580 //================================================================================
1582 * \brief Intersect two co-planar faces
1584 //================================================================================
1586 void Intersector::Algo::cutCoplanar()
1588 // find intersections of edges
1590 IntPoint2D intPoints[ 6 ];
1591 int nbIntPoints = 0;
1592 for ( int iE1 = 0; iE1 < 3; ++iE1 )
1594 int maxNbIntPoints = nbIntPoints + 2;
1595 for ( int iE2 = 0; iE2 < 3 && nbIntPoints < maxNbIntPoints; ++iE2 )
1596 nbIntPoints += intersectEdgeEdge( iE1, iE2, intPoints[ nbIntPoints ]);
1598 const int minNbOnPlane = Min( myNbOnPlane1, myNbOnPlane2 );
1600 if ( nbIntPoints == 0 ) // no intersections of edges
1603 if ( isPointInTriangle( myNodes1[0], myNodes2 )) // face2 includes face1
1605 else if ( isPointInTriangle( myNodes2[0], myNodes1 )) // face1 includes face2
1610 // add edges of an inner triangle to an outer one
1612 const std::vector< SMESH_NodeXYZ >& nodesIn = is1in2 ? myNodes1 : myNodes2;
1613 const SMDS_MeshElement* faceOut = is1in2 ? myFace2 : myFace1;
1614 const SMDS_MeshElement* faceIn = is1in2 ? myFace1 : myFace2;
1616 const CutFace& outFace = myCutFaces.Added( CutFace( faceOut ));
1617 CutLink link1( nodesIn.back().Node(), nodesIn.back().Node(), faceOut );
1618 CutLink link2( nodesIn.back().Node(), nodesIn.back().Node(), faceOut );
1620 link1.myIntNode = nodesIn.back();
1621 for ( size_t i = 0; i < nodesIn.size(); ++i )
1623 link2.myIntNode = nodesIn[ i ];
1624 outFace.AddEdge( link1, link2, faceIn, minNbOnPlane );
1625 link1.myIntNode = link2.myIntNode;
1630 // add parts of edges to a triangle including them
1632 CutLink link1, link2;
1633 IntPoint2D ip0, ip1;
1634 ip0.myU[0] = ip0.myU[1] = 0.;
1635 ip1.myU[0] = ip1.myU[1] = 1.;
1636 ip0.myEdgeInd[0] = ip0.myEdgeInd[1] = ip1.myEdgeInd[0] = ip1.myEdgeInd[1] = 0;
1638 for ( int isFromFace1 = 0; isFromFace1 < 2; ++isFromFace1 )
1640 const SMDS_MeshElement* faceTo = isFromFace1 ? myFace2 : myFace1;
1641 const SMDS_MeshElement* faceFrom = isFromFace1 ? myFace1 : myFace2;
1642 const std::vector< SMESH_NodeXYZ >& nodesTo = isFromFace1 ? myNodes2 : myNodes1;
1643 const std::vector< SMESH_NodeXYZ >& nodesFrom = isFromFace1 ? myNodes1 : myNodes2;
1644 const int iTo = isFromFace1 ? 1 : 0;
1645 const int iFrom = isFromFace1 ? 0 : 1;
1646 //const int nbOnPlaneFrom = isFromFace1 ? myNbOnPlane1 : myNbOnPlane2;
1648 const CutFace* cutFaceTo = & myCutFaces.Added( CutFace( faceTo ));
1649 // const CutFace* cutFaceFrom = 0;
1650 // if ( nbOnPlaneFrom > minNbOnPlane )
1651 // cutFaceFrom = & myCutFaces.Added( CutFace( faceTo ));
1653 link1.myFace = link2.myFace = faceTo;
1655 IntPoint2DCompare ipCompare( iFrom );
1656 TIntPointPtrSet pointsOnEdge( ipCompare ); // IntPoint2D sorted by parameter on edge
1658 for ( size_t iE = 0; iE < nodesFrom.size(); ++iE )
1660 // get parts of an edge iE
1662 ip0.myEdgeInd[ iTo ] = iE;
1663 ip1.myEdgeInd[ iTo ] = ( iE + 1 ) % nodesFrom.size();
1664 ip0.myNode = nodesFrom[ ip0.myEdgeInd[ iTo ]];
1665 ip1.myNode = nodesFrom[ ip1.myEdgeInd[ iTo ]];
1667 pointsOnEdge.clear();
1669 for ( int iP = 0; iP < nbIntPoints; ++iP )
1670 if ( intPoints[ iP ].myEdgeInd[ iFrom ] == iE )
1671 pointsOnEdge.insert( & intPoints[ iP ] );
1673 pointsOnEdge.insert( pointsOnEdge.begin(), & ip0 );
1674 pointsOnEdge.insert( pointsOnEdge.end(), & ip1 );
1676 // add edge parts to faceTo
1678 TIntPointPtrSet::iterator ipIt = pointsOnEdge.begin() + 1;
1679 for ( ; ipIt != pointsOnEdge.end(); ++ipIt )
1681 const IntPoint2D* p1 = *(ipIt-1);
1682 const IntPoint2D* p2 = *ipIt;
1683 gp_XYZ middle = 0.5 * ( p1->myNode + p2->myNode );
1684 if ( isPointInTriangle( middle, nodesTo ))
1686 p1->InitLink( link1, iTo, ( p1 != & ip0 ) ? nodesTo : nodesFrom );
1687 p2->InitLink( link2, iTo, ( p2 != & ip1 ) ? nodesTo : nodesFrom );
1688 cutFaceTo->AddEdge( link1, link2, faceFrom, minNbOnPlane );
1690 // if ( cutFaceFrom )
1692 // p1->InitLink( link1, iFrom, nodesFrom );
1693 // p2->InitLink( link2, iFrom, nodesFrom );
1694 // cutFaceTo->AddEdge( link1, link2, faceTo, minNbOnPlane );
1703 } // Intersector::Algo::cutCoplanar()
1705 //================================================================================
1707 * \brief Intersect edges added to myCutFaces
1709 //================================================================================
1711 void Intersector::Algo::IntersectNewEdges( const CutFace& cf )
1713 IntPoint2D intPoint;
1715 if ( cf.NbInternalEdges() < 2 )
1718 if ( myNodes1.empty() )
1724 const gp_XYZ& faceNorm = myNormals[ cf.myInitFace->GetID() ];
1725 setPlaneIndices( faceNorm ); // choose indices on an axis-aligned plane
1727 size_t limit = cf.myLinks.size() * cf.myLinks.size() * 2;
1730 while ( cf.myLinks[i1-1].IsInternal() && i1 > 0 )
1733 for ( ; i1 < cf.myLinks.size(); ++i1 )
1735 if ( !cf.myLinks[i1].IsInternal() )
1738 myIntPointSet.clear();
1739 for ( size_t i2 = i1 + 2; i2 < cf.myLinks.size(); ++i2 )
1741 if ( !cf.myLinks[i2].IsInternal() )
1744 // prepare to intersection
1745 myFace1 = cf.myLinks[i1].myFace;
1746 myNodes1[0] = cf.myLinks[i1].myNode1;
1747 myNodes1[1] = cf.myLinks[i1].myNode2;
1748 myFace2 = cf.myLinks[i2].myFace;
1749 myNodes2[0] = cf.myLinks[i2].myNode1;
1750 myNodes2[1] = cf.myLinks[i2].myNode2;
1753 intPoint.myIsCollinear = true; // to find collinear solutions
1754 if ( intersectEdgeEdge( 0, 0, intPoint ))
1756 if ( cf.myLinks[i1].IsSame( cf.myLinks[i2] )) // remove i2
1758 cf.myLinks[i1].ReplaceCoplanar( cf.myLinks[i2] );
1759 cf.myLinks.erase( cf.myLinks.begin() + i2, cf.myLinks.begin() + i2 + 2 );
1763 if ( !intPoint.myIsCollinear )
1765 intPoint.myEdgeInd[1] = i2;
1766 myIntPointSet.insert( intPoint );
1768 else // if ( intPoint.myIsCollinear ) // overlapping edges
1770 myIntPointSet.clear(); // to recompute
1772 if ( intPoint.myU[0] > intPoint.myU[1] ) // orient in same direction
1774 std::swap( intPoint.myU[0], intPoint.myU[1] );
1775 std::swap( myNodes1[0], myNodes1[1] );
1777 // replace _COPLANAR by _INTERNAL
1778 cf.myLinks[i1].ReplaceCoplanar( cf.myLinks[i1+1] );
1779 cf.myLinks[i2].ReplaceCoplanar( cf.myLinks[i2+1] );
1781 if ( coincide( myNodes1[0], myNodes2[0], myTol ) &&
1782 coincide( myNodes1[1], myNodes2[1], myTol ))
1784 cf.myLinks.erase( cf.myLinks.begin() + i2, cf.myLinks.begin() + i2 + 2 );
1789 EdgePart common = cf.myLinks[i1];
1790 common.ReplaceCoplanar( cf.myLinks[i2] );
1792 const SMDS_MeshNode* n1 = myNodes1[0].Node(); // end nodes of an overlapping part
1793 const SMDS_MeshNode* n2 = myNodes1[1].Node();
1794 size_t i3 = cf.myLinks.size();
1796 if ( myNodes1[0] != myNodes2[0] ) // a part before the overlapping one
1798 if ( intPoint.myU[0] < 0 )
1799 cf.myLinks[i1].Set( myNodes1[0].Node(), myNodes2[0].Node(),
1800 cf.myLinks[i1].myFace, cf.myLinks[i1].myIndex );
1802 cf.myLinks[i1].Set( myNodes2[0].Node(), myNodes1[0].Node(),
1803 cf.myLinks[i2].myFace, cf.myLinks[i2].myIndex );
1805 cf.myLinks[i1+1].Set( cf.myLinks[i1].myNode2,
1806 cf.myLinks[i1].myNode1,
1807 cf.myLinks[i1].myFace,
1808 cf.myLinks[i1].myIndex);
1809 n1 = cf.myLinks[i1].myNode2;
1814 if ( myNodes1[1] != myNodes2[1] ) // a part after the overlapping one
1816 if ( intPoint.myU[1] < 1 )
1817 cf.myLinks[i2].Set( myNodes1[1].Node(), myNodes2[1].Node(),
1818 cf.myLinks[i2].myFace, cf.myLinks[i2].myIndex );
1820 cf.myLinks[i2].Set( myNodes2[1].Node(), myNodes1[1].Node(),
1821 cf.myLinks[i1].myFace, cf.myLinks[i1].myIndex );
1823 cf.myLinks[i2+1].Set( cf.myLinks[i2].myNode2,
1824 cf.myLinks[i2].myNode1,
1825 cf.myLinks[i2].myFace,
1826 cf.myLinks[i2].myIndex);
1827 n2 = cf.myLinks[i2].myNode1;
1832 if ( i3 == cf.myLinks.size() )
1833 cf.myLinks.resize( i3 + 2 );
1835 cf.myLinks[i3].Set ( n1, n2, common.myFace, common.myIndex );
1836 cf.myLinks[i3+1].Set( n2, n1, common.myFace, common.myIndex );
1838 i2 = i1 + 1; // recheck modified i1
1843 // // remember a new node
1844 // CutLink link1( myNodes1[0].Node(), myNodes1[1].Node(), cf.myInitFace );
1845 // CutLink link2( myNodes2[0].Node(), myNodes2[1].Node(), cf.myInitFace );
1846 // link2.myIntNode = link1.myIntNode = intPoint.myNode;
1847 // addLink( link1 );
1848 // addLink( link2 );
1851 // size_t i = cf.myLinks.size();
1852 // if ( intPoint.myNode != cf.myLinks[ i1 ].myNode1 &&
1853 // intPoint.myNode != cf.myLinks[ i1 ].myNode2 )
1855 // cf.myLinks.push_back( cf.myLinks[ i1 ]);
1856 // cf.myLinks.push_back( cf.myLinks[ i1 + 1 ]);
1857 // cf.myLinks[ i1 ].myNode2 = cf.myLinks[ i1 + 1 ].myNode1 = intPoint.Node();
1858 // cf.myLinks[ i ].myNode1 = cf.myLinks[ i + 1 ].myNode2 = intPoint.Node();
1860 // if ( intPoint.myNode != cf.myLinks[ i2 ].myNode1 &&
1861 // intPoint.myNode != cf.myLinks[ i2 ].myNode2 )
1863 // i = cf.myLinks.size();
1864 // cf.myLinks.push_back( cf.myLinks[ i2 ]);
1865 // cf.myLinks.push_back( cf.myLinks[ i2 + 1 ]);
1866 // cf.myLinks[ i2 ].myNode2 = cf.myLinks[ i2 + 1 ].myNode1 = intPoint.Node();
1867 // cf.myLinks[ i ].myNode1 = cf.myLinks[ i + 1 ].myNode2 = intPoint.Node();
1871 } // if ( intersectEdgeEdge( 0, 0, intPoint ))
1877 // split i1 edge and all edges it intersects
1878 // don't do it inside intersection loop in order not to loose direction of i1 edge
1879 if ( !myIntPointSet.empty() )
1881 cf.myLinks.reserve( cf.myLinks.size() + myIntPointSet.size() * 2 + 2 );
1883 EdgePart* edge1 = &cf.myLinks[ i1 ];
1884 EdgePart* twin1 = &cf.myLinks[ i1 + 1 ];
1886 TIntPointSet::iterator ipIt = myIntPointSet.begin();
1887 for ( ; ipIt != myIntPointSet.end(); ++ipIt ) // int points sorted on i1 edge
1889 size_t i = cf.myLinks.size();
1890 if ( ipIt->myNode != edge1->myNode1 &&
1891 ipIt->myNode != edge1->myNode2 )
1893 cf.myLinks.push_back( *edge1 );
1894 cf.myLinks.push_back( *twin1 );
1895 edge1->myNode2 = twin1->myNode1 = ipIt->Node();
1896 cf.myLinks[ i ].myNode1 = cf.myLinks[ i + 1 ].myNode2 = ipIt->Node();
1897 edge1 = & cf.myLinks[ i ];
1898 twin1 = & cf.myLinks[ i + 1 ];
1900 size_t i2 = ipIt->myEdgeInd[1];
1901 if ( ipIt->myNode != cf.myLinks[ i2 ].myNode1 &&
1902 ipIt->myNode != cf.myLinks[ i2 ].myNode2 )
1904 i = cf.myLinks.size();
1905 cf.myLinks.push_back( cf.myLinks[ i2 ]);
1906 cf.myLinks.push_back( cf.myLinks[ i2 + 1 ]);
1907 cf.myLinks[ i2 ].myNode2 = cf.myLinks[ i2 + 1 ].myNode1 = ipIt->Node();
1908 cf.myLinks[ i ].myNode1 = cf.myLinks[ i + 1 ].myNode2 = ipIt->Node();
1911 if ( cf.myLinks.size() >= limit )
1912 throw SALOME_Exception( "Infinite loop in Intersector::Algo::IntersectNewEdges()" );
1914 ++i1; // each internal edge encounters twice
1919 //================================================================================
1921 * \brief Split intersected faces
1923 //================================================================================
1925 void Intersector::Algo::MakeNewFaces( SMESH_MeshAlgos::TElemIntPairVec& theNew2OldFaces,
1926 SMESH_MeshAlgos::TNodeIntPairVec& theNew2OldNodes,
1927 const double theSign,
1928 const bool theOptimize)
1930 // fill theNew2OldFaces if empty
1931 TCutFaceMap::const_iterator cutFacesIt = myCutFaces.cbegin();
1932 if ( theNew2OldFaces.empty() )
1933 for ( ; cutFacesIt != myCutFaces.cend(); ++cutFacesIt )
1935 const CutFace& cf = *cutFacesIt;
1936 int index = cf.myInitFace->GetID(); // index in theNew2OldFaces
1937 if ((int) theNew2OldFaces.size() <= index )
1938 theNew2OldFaces.resize( index + 1 );
1939 theNew2OldFaces[ index ] = std::make_pair( cf.myInitFace, index );
1942 // unmark all nodes except intersection ones
1944 for ( SMDS_NodeIteratorPtr nIt = myMesh->nodesIterator(); nIt->more(); )
1946 const SMDS_MeshNode* n = nIt->next();
1947 if ( n->isMarked() && n->GetID()-1 < (int) theNew2OldNodes.size() )
1948 n->setIsMarked( false );
1950 // SMESH_MeshAlgos::MarkElems( myMesh->nodesIterator(), false );
1952 TCutLinkMap::const_iterator cutLinksIt = myCutLinks.cbegin();
1953 // for ( ; cutLinksIt != myCutLinks.cend(); ++cutLinksIt )
1955 // const CutLink& link = *cutLinksIt;
1956 // if ( link.IntNode() && link.IntNode()->GetID()-1 < (int) theNew2OldNodes.size() )
1957 // link.IntNode()->setIsMarked( true );
1960 // intersect edges added to myCutFaces
1962 for ( cutFacesIt = myCutFaces.cbegin(); cutFacesIt != myCutFaces.cend(); ++cutFacesIt )
1964 const CutFace& cf = *cutFacesIt;
1965 cf.ReplaceNodes( myRemove2KeepNodes );
1966 IntersectNewEdges( cf );
1971 EdgeLoopSet loopSet;
1972 SMESH_MeshAlgos::Triangulate triangulator( theOptimize );
1973 std::vector< EdgePart > cutOffLinks;
1974 TLinkMap cutOffCoplanarLinks;
1975 std::vector< const CutFace* > touchedFaces;
1976 SMESH_MeshAlgos::TElemIntPairVec::value_type new2OldTria;
1978 std::vector< const SMDS_MeshNode* > nodes;
1979 std::vector<const SMDS_MeshElement *> faces;
1981 cutOffLinks.reserve( myCutFaces.Extent() * 2 );
1983 for ( cutFacesIt = myCutFaces.cbegin(); cutFacesIt != myCutFaces.cend(); ++cutFacesIt )
1985 const CutFace& cf = *cutFacesIt;
1988 touchedFaces.push_back( & cf );
1992 const gp_XYZ& normal = myNormals[ cf.myInitFace->GetID() ];
1994 // form loops of new faces
1995 cf.ReplaceNodes( myRemove2KeepNodes );
1996 cf.MakeLoops( loopSet, normal );
1998 // avoid loops that are not connected to boundary edges of cf.myInitFace
1999 if ( cf.RemoveInternalLoops( loopSet ))
2001 IntersectNewEdges( cf );
2002 cf.MakeLoops( loopSet, normal );
2004 // erase loops that are cut off by face intersections
2005 cf.CutOffLoops( loopSet, theSign, myNormals, cutOffLinks, cutOffCoplanarLinks );
2007 int index = cf.myInitFace->GetID(); // index in theNew2OldFaces
2009 const SMDS_MeshElement* tria;
2010 for ( size_t iL = 0; iL < loopSet.myNbLoops; ++iL )
2012 EdgeLoop& loop = loopSet.myLoops[ iL ];
2013 if ( loop.myLinks.size() == 0 )
2016 int nbTria = triangulator.GetTriangles( &loop, nodes );
2017 int nbNodes = 3 * nbTria;
2018 for ( int i = 0; i < nbNodes; i += 3 )
2020 if ( nodes[i] == nodes[i+1] || nodes[i] == nodes[i+2] || nodes[i+1] == nodes[i+2] )
2023 std::cerr << "BAD tria" << std::endl;
2028 if (!( tria = myMesh->FindFace( nodes[i], nodes[i+1], nodes[i+2] )))
2029 tria = myMesh->AddFace( nodes[i], nodes[i+1], nodes[i+2] );
2030 tria->setIsMarked( true ); // not to remove it
2032 new2OldTria = std::make_pair( tria, theNew2OldFaces[ index ].second );
2033 if ( tria->GetID() < (int)theNew2OldFaces.size() )
2034 theNew2OldFaces[ tria->GetID() ] = new2OldTria;
2036 theNew2OldFaces.push_back( new2OldTria );
2038 if ( index == tria->GetID() )
2039 index = 0; // do not remove tria
2042 theNew2OldFaces[ index ].first = 0;
2045 // remove split faces
2046 for ( size_t id = 1; id < theNew2OldFaces.size(); ++id )
2048 if ( theNew2OldFaces[id].first ||
2049 theNew2OldFaces[id].second == 0 )
2051 if ( const SMDS_MeshElement* f = myMesh->FindElement( id ))
2052 myMesh->RemoveFreeElement( f );
2055 // remove faces that are merged off
2056 for ( cutFacesIt = myCutFaces.cbegin(); cutFacesIt != myCutFaces.cend(); ++cutFacesIt )
2058 const CutFace& cf = *cutFacesIt;
2059 if ( !cf.myLinks.empty() || cf.myInitFace->IsNull() )
2062 nodes.assign( cf.myInitFace->begin_nodes(), cf.myInitFace->end_nodes() );
2063 for ( size_t i = 0; i < nodes.size(); ++i )
2065 const SMDS_MeshNode* n = nodes[ i ];
2066 while ( myRemove2KeepNodes.IsBound( n ))
2067 n = myRemove2KeepNodes( n );
2068 if ( n != nodes[ i ] && cf.myInitFace->GetNodeIndex( n ) >= 0 )
2070 theNew2OldFaces[ cf.myInitFace->GetID() ].first = 0;
2071 myMesh->RemoveFreeElement( cf.myInitFace );
2077 // remove faces connected to cut off parts of cf.myInitFace
2080 for ( size_t i = 0; i < cutOffLinks.size(); ++i )
2083 nodes[0] = cutOffLinks[i].myNode1;
2084 nodes[1] = cutOffLinks[i].myNode2;
2086 if ( nodes[0] != nodes[1] &&
2087 myMesh->GetElementsByNodes( nodes, faces ))
2089 if ( // cutOffLinks[i].myFace &&
2090 cutOffLinks[i].myIndex != EdgePart::_COPLANAR &&
2093 for ( size_t iF = 0; iF < faces.size(); ++iF )
2095 int index = faces[iF]->GetID();
2096 // if ( //faces[iF]->isMarked() || // kept part of cutFace
2097 // !theNew2OldFaces[ index ].first ) // already removed
2099 cutFace.myInitFace = faces[iF];
2100 // if ( myCutFaces.Contains( cutFace )) // keep cutting faces needed in CutOffLoops()
2102 // if ( !myCutFaces.Added( cutFace ).IsCut() )
2103 // theNew2OldFaces[ index ].first = 0;
2106 cutFace.myLinks.clear();
2107 cutFace.InitLinks();
2108 for ( size_t iL = 0; iL < cutFace.myLinks.size(); ++iL )
2109 if ( !cutOffLinks[i].IsSame( cutFace.myLinks[ iL ]))
2110 cutOffLinks.push_back( cutFace.myLinks[ iL ]);
2112 theNew2OldFaces[ index ].first = 0;
2113 myMesh->RemoveFreeElement( faces[iF] );
2118 // replace nodes in touched faces
2120 // treat touched faces
2121 for ( size_t i = 0; i < touchedFaces.size(); ++i )
2123 const CutFace& cf = *touchedFaces[i];
2124 if ( cf.myInitFace->IsNull() )
2127 int index = cf.myInitFace->GetID(); // index in theNew2OldFaces
2128 if ( !theNew2OldFaces[ index ].first )
2129 continue; // already cut off
2132 if ( !cf.ReplaceNodes( myRemove2KeepNodes ))
2134 if ( cf.myLinks.size() == 3 &&
2135 cf.myInitFace->GetNodeIndex( cf.myLinks[0].myNode1 ) >= 0 &&
2136 cf.myInitFace->GetNodeIndex( cf.myLinks[1].myNode1 ) >= 0 &&
2137 cf.myInitFace->GetNodeIndex( cf.myLinks[2].myNode1 ) >= 0 )
2138 continue; // just keep as is
2141 if ( cf.myLinks.size() == 3 )
2143 const SMDS_MeshElement* tria = myMesh->AddFace( cf.myLinks[0].myNode1,
2144 cf.myLinks[1].myNode1,
2145 cf.myLinks[2].myNode1 );
2146 new2OldTria = std::make_pair( tria, theNew2OldFaces[ index ].second );
2147 if ( tria->GetID() < (int)theNew2OldFaces.size() )
2148 theNew2OldFaces[ tria->GetID() ] = new2OldTria;
2150 theNew2OldFaces.push_back( new2OldTria );
2152 theNew2OldFaces[ index ].first = 0;
2156 // add used new nodes to theNew2OldNodes
2157 SMESH_MeshAlgos::TNodeIntPairVec::value_type new2OldNode;
2158 new2OldNode.second = 0;
2159 for ( cutLinksIt = myCutLinks.cbegin(); cutLinksIt != myCutLinks.cend(); ++cutLinksIt )
2161 const CutLink& link = *cutLinksIt;
2162 if ( link.IntNode() ) // && link.IntNode()->NbInverseElements() > 0 )
2164 new2OldNode.first = link.IntNode();
2165 theNew2OldNodes.push_back( new2OldNode );
2172 //================================================================================
2173 Intersector::Intersector( SMDS_Mesh* mesh, double tol, const std::vector< gp_XYZ >& normals )
2175 myAlgo = new Algo( mesh, tol, normals );
2177 //================================================================================
2178 Intersector::~Intersector()
2182 //================================================================================
2183 //! compute cut of two faces of the mesh
2184 void Intersector::Cut( const SMDS_MeshElement* face1,
2185 const SMDS_MeshElement* face2,
2186 const int nbCommonNodes )
2188 myAlgo->Cut( face1, face2, nbCommonNodes );
2190 //================================================================================
2191 //! store a face cut by a line given by its ends
2192 // accompanied by indices of intersected face edges.
2193 // Edge index is <0 if a line end is inside the face.
2194 void Intersector::Cut( const SMDS_MeshElement* face,
2195 SMESH_NodeXYZ& lineEnd1,
2197 SMESH_NodeXYZ& lineEnd2,
2200 myAlgo->Cut( face, lineEnd1, edgeIndex1, lineEnd2, edgeIndex2 );
2202 //================================================================================
2203 //! split all face intersected by Cut() methods
2204 void Intersector::MakeNewFaces( SMESH_MeshAlgos::TElemIntPairVec& theNew2OldFaces,
2205 SMESH_MeshAlgos::TNodeIntPairVec& theNew2OldNodes,
2206 const double theSign,
2207 const bool theOptimize )
2209 myAlgo->MakeNewFaces( theNew2OldFaces, theNew2OldNodes, theSign, theOptimize );
2211 //================================================================================
2212 //! Cut a face by planes, whose normals point to parts to keep
2213 bool Intersector::CutByPlanes(const SMDS_MeshElement* theFace,
2214 const std::vector< gp_Ax1 > & thePlanes,
2215 const double theTol,
2216 std::vector< TFace > & theNewFaceConnectivity )
2218 theNewFaceConnectivity.clear();
2220 // check if theFace is wholly cut off
2221 std::vector< SMESH_NodeXYZ > facePoints( theFace->begin_nodes(), theFace->end_nodes() );
2222 facePoints.resize( theFace->NbCornerNodes() );
2223 for ( size_t iP = 0; iP < thePlanes.size(); ++iP )
2226 const gp_Pnt& O = thePlanes[iP].Location();
2227 for ( size_t i = 0; i < facePoints.size(); ++i )
2229 gp_Vec Op( O, facePoints[i] );
2230 nbOut += ( Op * thePlanes[iP].Direction() <= 0 );
2232 if ( nbOut == facePoints.size() )
2236 // copy theFace into a temporary mesh
2239 std::vector< const SMDS_MeshNode* > faceNodes;
2240 faceNodes.resize( facePoints.size() );
2241 for ( size_t i = 0; i < facePoints.size(); ++i )
2243 const SMESH_NodeXYZ& n = facePoints[i];
2244 faceNodes[i] = mesh.AddNode( n.X(), n.Y(), n.Z() );
2247 const SMDS_MeshElement* faceToCut = 0;
2248 switch ( theFace->NbCornerNodes() )
2251 faceToCut = mesh.AddFace( faceNodes[0], faceNodes[1], faceNodes[2] );
2254 faceToCut = mesh.AddFace( faceNodes[0], faceNodes[1], faceNodes[2], faceNodes[3] );
2257 faceToCut = mesh.AddPolygonalFace( faceNodes );
2260 std::vector< gp_XYZ > normals( 2 + thePlanes.size() );
2261 SMESH_MeshAlgos::FaceNormal( faceToCut, normals[ faceToCut->GetID() ]);
2263 // add faces corresponding to thePlanes
2264 std::vector< const SMDS_MeshElement* > planeFaces;
2265 double faceSize = Sqrt( faceBox.SquareExtent() );
2266 gp_XYZ center = 0.5 * ( faceBox.CornerMin() + faceBox.CornerMax() );
2267 for ( size_t i = 0; i < thePlanes.size(); ++i )
2269 gp_Ax2 plnAx( thePlanes[i].Location(), thePlanes[i].Direction() );
2270 gp_XYZ O = plnAx.Location().XYZ();
2271 gp_XYZ X = plnAx.XDirection().XYZ();
2272 gp_XYZ Y = plnAx.YDirection().XYZ();
2273 gp_XYZ Z = plnAx.Direction().XYZ();
2275 double dot = ( O - center ) * Z;
2276 gp_XYZ o = center + Z * dot; // center projected to a plane
2278 gp_XYZ p1 = o + X * faceSize * 2;
2279 gp_XYZ p2 = o + Y * faceSize * 2;
2280 gp_XYZ p3 = o - (X + Y ) * faceSize * 2;
2282 const SMDS_MeshNode* n1 = mesh.AddNode( p1.X(), p1.Y(), p1.Z() );
2283 const SMDS_MeshNode* n2 = mesh.AddNode( p2.X(), p2.Y(), p2.Z() );
2284 const SMDS_MeshNode* n3 = mesh.AddNode( p3.X(), p3.Y(), p3.Z() );
2285 planeFaces.push_back( mesh.AddFace( n1, n2, n3 ));
2287 normals[ planeFaces.back()->GetID() ] = thePlanes[i].Direction().XYZ();
2291 Algo algo ( &mesh, theTol, normals );
2292 for ( size_t i = 0; i < planeFaces.size(); ++i )
2294 algo.Cut( faceToCut, planeFaces[i], 0 );
2297 // retrieve a result
2298 SMESH_MeshAlgos::TElemIntPairVec new2OldFaces;
2299 SMESH_MeshAlgos::TNodeIntPairVec new2OldNodes;
2300 TCutFaceMap::const_iterator cutFacesIt= algo.myCutFaces.cbegin();
2301 for ( ; cutFacesIt != algo.myCutFaces.cend(); ++cutFacesIt )
2303 const CutFace& cf = *cutFacesIt;
2304 if ( cf.myInitFace != faceToCut )
2309 theNewFaceConnectivity.push_back( facePoints );
2313 // intersect cut lines
2314 algo.IntersectNewEdges( cf );
2316 // form loops of new faces
2317 EdgeLoopSet loopSet;
2318 cf.MakeLoops( loopSet, normals[ faceToCut->GetID() ]);
2320 // erase loops that are cut off by thePlanes
2321 const double sign = 1;
2322 std::vector< EdgePart > cutOffLinks;
2323 TLinkMap cutOffCoplanarLinks;
2324 cf.CutOffLoops( loopSet, sign, normals, cutOffLinks, cutOffCoplanarLinks );
2326 for ( size_t iL = 0; iL < loopSet.myNbLoops; ++iL )
2328 EdgeLoop& loop = loopSet.myLoops[ iL ];
2329 if ( loop.myLinks.size() > 0 )
2332 for ( SMDS_NodeIteratorPtr nIt = loop.nodeIterator(); nIt->more(); )
2334 const SMDS_MeshNode* n = nIt->next();
2335 facePoints.push_back( n );
2336 int iN = faceToCut->GetNodeIndex( n );
2338 facePoints.back()._node = 0; // an intersection point
2340 facePoints.back()._node = theFace->GetNode( iN );
2342 theNewFaceConnectivity.push_back( facePoints );
2348 return theNewFaceConnectivity.empty();
2351 } // namespace SMESH_MeshAlgos
2355 //================================================================================
2359 //================================================================================
2361 void CutFace::Dump() const
2363 std::cout << std::endl << "INI F " << myInitFace->GetID() << std::endl;
2364 for ( size_t i = 0; i < myLinks.size(); ++i )
2365 std::cout << "[" << i << "] ("
2366 << char(( myLinks[i].IsInternal() ? 'j' : '0' ) + myLinks[i].myIndex ) << ") "
2367 << myLinks[i].myNode1->GetID() << " - " << myLinks[i].myNode2->GetID()
2368 << " " << ( myLinks[i].myFace ? 'F' : 'C' )
2369 << ( myLinks[i].myFace ? myLinks[i].myFace->GetID() : 0 ) << " " << std::endl;
2372 //================================================================================
2374 * \brief Add an edge cutting this face
2375 * \param [in] p1 - start point of the edge
2376 * \param [in] p2 - end point of the edge
2377 * \param [in] cutter - a face producing the added cut edge.
2378 * \param [in] nbOnPlane - nb of triangle nodes lying on the plane of the cutter face
2380 //================================================================================
2382 void CutFace::AddEdge( const CutLink& p1,
2384 const SMDS_MeshElement* cutterFace,
2385 const int nbOnPlane,
2386 const int iNotOnPlane) const
2388 int iN[2] = { myInitFace->GetNodeIndex( p1.IntNode() ),
2389 myInitFace->GetNodeIndex( p2.IntNode() ) };
2390 if ( iN[0] >= 0 && iN[1] >= 0 )
2392 // the cutting edge is a whole side
2393 if (( cutterFace && nbOnPlane < 3 ) &&
2394 !( cutterFace->GetNodeIndex( p1.IntNode() ) >= 0 &&
2395 cutterFace->GetNodeIndex( p2.IntNode() ) >= 0 ))
2398 myLinks[ Abs( iN[0] - iN[1] ) == 1 ? Min( iN[0], iN[1] ) : 2 ].myFace = cutterFace;
2403 if ( p1.IntNode() == p2.IntNode() )
2405 AddPoint( p1, p2, 1e-10 );
2411 // cut side edges by a new one
2413 int iEOnPlane = ( nbOnPlane == 2 ) ? ( iNotOnPlane + 1 ) % 3 : -1;
2416 for ( int is2nd = 0; is2nd < 2; ++is2nd )
2418 const CutLink& p = is2nd ? p2 : p1;
2420 if ( iN[ is2nd ] >= 0 )
2423 int iE = Max( iEOnPlane, myInitFace->GetNodeIndex( p.Node1() ));
2425 continue; // link of other face
2427 SMESH_NodeXYZ n0 = myLinks[iE].myNode1;
2428 dist[ is2nd ] = ( n0 - p.myIntNode ).SquareModulus();
2430 for ( size_t i = 0; i < myLinks.size(); ++i )
2431 if ( myLinks[i].myIndex == iE )
2433 double d1 = n0.SquareDistance( myLinks[i].myNode1 );
2434 if ( d1 < dist[ is2nd ] )
2436 double d2 = n0.SquareDistance( myLinks[i].myNode2 );
2437 if ( dist[ is2nd ] < d2 )
2439 myLinks.push_back( myLinks[i] );
2440 myLinks.back().myNode1 = myLinks[i].myNode2 = p.IntNode();
2447 int state = nbOnPlane == 3 ? EdgePart::_COPLANAR : EdgePart::_INTERNAL;
2449 // look for an existing equal edge
2450 if ( nbOnPlane == 2 )
2452 SMESH_NodeXYZ n0 = myLinks[ iEOnPlane ].myNode1;
2453 if ( iN[0] >= 0 ) dist[0] = ( n0 - p1.myIntNode ).SquareModulus();
2454 if ( iN[1] >= 0 ) dist[1] = ( n0 - p2.myIntNode ).SquareModulus();
2455 if ( dist[0] > dist[1] )
2456 std::swap( dist[0], dist[1] );
2457 for ( size_t i = 0; i < myLinks.size(); ++i )
2459 if ( myLinks[i].myIndex != iEOnPlane )
2461 gp_XYZ mid = 0.5 * ( SMESH_NodeXYZ( myLinks[i].myNode1 ) +
2462 SMESH_NodeXYZ( myLinks[i].myNode2 ));
2463 double d = ( n0 - mid ).SquareModulus();
2464 if ( dist[0] < d && d < dist[1] )
2465 myLinks[i].myFace = cutterFace;
2471 EdgePart newEdge; newEdge.Set( p1.IntNode(), p2.IntNode(), cutterFace, state );
2472 for ( size_t i = 0; i < myLinks.size(); ++i )
2474 if ( myLinks[i].IsSame( newEdge ))
2476 // if ( !myLinks[i].IsInternal() )
2477 // myLinks[ i ].myFace = cutterFace;
2479 myLinks[ i ].ReplaceCoplanar( newEdge );
2480 if ( myLinks[i].IsInternal() && i+1 < myLinks.size() )
2481 myLinks[ i+1 ].ReplaceCoplanar( newEdge );
2484 i += myLinks[i].IsInternal();
2488 size_t i = myLinks.size();
2489 myLinks.resize( i + 2 );
2490 myLinks[ i ].Set( p1.IntNode(), p2.IntNode(), cutterFace, state );
2491 myLinks[ i+1 ].Set( p2.IntNode(), p1.IntNode(), cutterFace, state );
2494 //================================================================================
2496 * \brief Add a point cutting this face
2498 //================================================================================
2500 void CutFace::AddPoint( const CutLink& p1, const CutLink& p2, double tol ) const
2502 if ( myInitFace->GetNodeIndex( p1.IntNode() ) >= 0 ||
2503 myInitFace->GetNodeIndex( p2.IntNode() ) >= 0 )
2508 const CutLink* link = &p1;
2509 int iE = myInitFace->GetNodeIndex( link->Node1() );
2513 iE = myInitFace->GetNodeIndex( link->Node1() );
2517 // cut an existing edge by the point
2518 SMESH_NodeXYZ n0 = link->Node1();
2519 double d = ( n0 - link->myIntNode ).SquareModulus();
2521 for ( size_t i = 0; i < myLinks.size(); ++i )
2522 if ( myLinks[i].myIndex == iE )
2524 double d1 = n0.SquareDistance( myLinks[i].myNode1 );
2527 double d2 = n0.SquareDistance( myLinks[i].myNode2 );
2530 myLinks.push_back( myLinks[i] );
2531 myLinks.back().myNode1 = myLinks[i].myNode2 = link->IntNode();
2537 else // point is inside the triangle
2539 // // check if a point already added
2540 // for ( size_t i = 3; i < myLinks.size(); ++i )
2541 // if ( myLinks[i].myNode1 == p1.IntNode() )
2544 // // create a link between the point and the closest corner node
2545 // const SMDS_MeshNode* closeNode = myLinks[0].myNode1;
2546 // double minDist = p1.myIntNode.SquareDistance( closeNode );
2547 // for ( int i = 1; i < 3; ++i )
2549 // double dist = p1.myIntNode.SquareDistance( myLinks[i].myNode1 );
2550 // if ( dist < minDist )
2553 // closeNode = myLinks[i].myNode1;
2556 // if ( minDist > tol * tol )
2558 // size_t i = myLinks.size();
2559 // myLinks.resize( i + 2 );
2560 // myLinks[ i ].Set( p1.IntNode(), closeNode );
2561 // myLinks[ i+1 ].Set( closeNode, p1.IntNode() );
2566 //================================================================================
2568 * \brief Perform node replacement
2570 //================================================================================
2572 bool CutFace::ReplaceNodes( const TNNMap& theRm2KeepMap ) const
2574 bool replaced = false;
2575 for ( size_t i = 0; i < myLinks.size(); ++i )
2577 while ( theRm2KeepMap.IsBound( myLinks[i].myNode1 ))
2578 replaced = ( myLinks[i].myNode1 = theRm2KeepMap( myLinks[i].myNode1 ));
2580 while ( theRm2KeepMap.IsBound( myLinks[i].myNode2 ))
2581 replaced = ( myLinks[i].myNode2 = theRm2KeepMap( myLinks[i].myNode2 ));
2584 //if ( replaced ) // remove equal links
2586 for ( size_t i1 = 0; i1 < myLinks.size(); ++i1 )
2588 if ( myLinks[i1].myNode1 == myLinks[i1].myNode2 )
2590 myLinks.erase( myLinks.begin() + i1,
2591 myLinks.begin() + i1 + 1 + myLinks[i1].IsInternal() );
2595 size_t i2 = i1 + 1 + myLinks[i1].IsInternal();
2596 for ( ; i2 < myLinks.size(); ++i2 )
2598 if ( !myLinks[i2].IsInternal() )
2600 if ( myLinks[i1].IsSame( myLinks[i2] ))
2602 myLinks[i1]. ReplaceCoplanar( myLinks[i2] );
2603 if ( myLinks[i1].IsInternal() )
2604 myLinks[i1+1].ReplaceCoplanar( myLinks[i2+1] );
2605 if ( !myLinks[i1].myFace && myLinks[i2].myFace )
2607 myLinks[i1]. myFace = myLinks[i2].myFace;
2608 if ( myLinks[i1].IsInternal() )
2609 myLinks[i1+1].myFace = myLinks[i2+1].myFace;
2611 myLinks.erase( myLinks.begin() + i2,
2612 myLinks.begin() + i2 + 2 );
2618 i1 += myLinks[i1].IsInternal();
2625 //================================================================================
2627 * \brief Initialize myLinks with edges of myInitFace
2629 //================================================================================
2631 void CutFace::InitLinks() const
2633 if ( !myLinks.empty() ) return;
2635 int nbNodes = myInitFace->NbNodes();
2636 myLinks.reserve( nbNodes * 2 );
2637 myLinks.resize( nbNodes );
2639 for ( int i = 0; i < nbNodes; ++i )
2641 const SMDS_MeshNode* n1 = myInitFace->GetNode( i );
2642 const SMDS_MeshNode* n2 = myInitFace->GetNodeWrap( i + 1);
2643 myLinks[i].Set( n1, n2, 0, i );
2647 //================================================================================
2649 * \brief Return number of internal edges
2651 //================================================================================
2653 int CutFace::NbInternalEdges() const
2656 for ( size_t i = 3; i < myLinks.size(); ++i )
2657 nb += myLinks[i].IsInternal();
2659 return nb / 2; // each internal edge encounters twice
2662 //================================================================================
2664 * \brief Remove loops that are not connected to boundary edges of myFace by
2665 * adding edges connecting these loops to the boundary.
2666 * Such loops must be removed as they form polygons with ring topology.
2668 //================================================================================
2670 bool CutFace::RemoveInternalLoops( EdgeLoopSet& theLoops ) const
2672 size_t nbReachedLoops = 0;
2674 // count loops including boundary EdgeParts
2675 for ( size_t iL = 0; iL < theLoops.myNbLoops; ++iL )
2677 EdgeLoop& loop = theLoops.myLoops[ iL ];
2679 for ( size_t iE = 0; iE < loop.myLinks.size(); ++iE )
2680 if ( !loop.myLinks[ iE ]->IsInternal() )
2682 nbReachedLoops += loop.SetConnected();
2686 if ( nbReachedLoops == theLoops.myNbLoops )
2687 return false; // no unreachable loops
2690 // try to reach all loops by propagating via internal edges shared by loops
2691 size_t prevNbReached;
2694 prevNbReached = nbReachedLoops;
2696 for ( size_t iL = 0; iL < theLoops.myNbLoops; ++iL )
2698 EdgeLoop& loop = theLoops.myLoops[ iL ];
2699 if ( !loop.myIsBndConnected )
2702 for ( size_t iE = 0; iE < loop.myLinks.size(); ++iE )
2703 if ( loop.myLinks[ iE ]->IsInternal() )
2705 const EdgePart* twinEdge = getTwin( loop.myLinks[ iE ]);
2706 EdgeLoop* loop2 = theLoops.GetLoopOf( twinEdge );
2707 if ( loop2->SetConnected() && ++nbReachedLoops == theLoops.myNbLoops )
2708 return false; // no unreachable loops
2712 while ( prevNbReached < nbReachedLoops );
2716 for ( size_t iL = 0; iL < theLoops.myNbLoops; ++iL )
2718 EdgeLoop& loop = theLoops.myLoops[ iL ];
2719 if ( loop.myIsBndConnected || loop.myLinks.size() == 0 )
2722 if ( loop.myHasPending )
2724 // try to join the loop to another one, with which it contacts at a node
2726 // look for a node where the loop reverses
2727 const EdgePart* edgePrev = loop.myLinks.back();
2728 for ( size_t iE = 0; iE < loop.myLinks.size(); edgePrev = loop.myLinks[ iE++ ] )
2730 if ( !edgePrev->IsTwin( *loop.myLinks[ iE ]))
2732 const SMDS_MeshNode* reverseNode = edgePrev->myNode2;
2734 // look for a loop including reverseNode
2735 size_t iContactEdge2; // index(+1) of edge starting at reverseNode
2736 for ( size_t iL2 = 0; iL2 < theLoops.myNbLoops; ++iL2 )
2740 EdgeLoop& loop2 = theLoops.myLoops[ iL2 ];
2741 if ( ! ( iContactEdge2 = loop2.Contains( reverseNode )))
2744 // insert loop2 into the loop
2745 theLoops.Join( loop, iE, loop2, iContactEdge2 - 1 );
2748 if ( loop.myIsBndConnected )
2752 if ( loop.myIsBndConnected )
2756 // add links connecting internal loops with the boundary ones
2758 // find a pair of closest nodes
2759 const SMDS_MeshNode *closestNode1, *closestNode2;
2760 double minDist = 1e100;
2761 for ( size_t iE = 0; iE < loop.myLinks.size(); ++iE )
2763 SMESH_NodeXYZ n1 = loop.myLinks[ iE ]->myNode1;
2765 for ( size_t i = 0; i < myLinks.size(); ++i )
2767 if ( !loop.Contains( myLinks[i].myNode1 ))
2769 double dist = n1.SquareDistance( myLinks[i].myNode1 );
2770 if ( dist < minDist )
2773 closestNode1 = loop.myLinks[ iE ]->myNode1;
2774 closestNode2 = myLinks[i].myNode1;
2777 if ( myLinks[i].IsInternal() )
2782 size_t i = myLinks.size();
2783 myLinks.resize( i + 2 );
2784 myLinks[ i ].Set( closestNode1, closestNode2 );
2785 myLinks[ i+1 ].Set( closestNode2, closestNode1 );
2791 //================================================================================
2793 * \brief Return equal reversed edge
2795 //================================================================================
2797 EdgePart* CutFace::getTwin( const EdgePart* edge ) const
2799 size_t i = edge - & myLinks[0];
2801 if ( i > 2 && myLinks[ i-1 ].IsTwin( *edge ))
2802 return & myLinks[ i-1 ];
2804 if ( i+1 < myLinks.size() &&
2805 myLinks[ i+1 ].IsTwin( *edge ))
2806 return & myLinks[ i+1 ];
2811 //================================================================================
2813 * \brief Fill loops of edges
2815 //================================================================================
2817 void CutFace::MakeLoops( EdgeLoopSet& theLoops, const gp_XYZ& theFaceNorm ) const
2819 theLoops.Init( myLinks );
2821 if ( myLinks.size() == 3 )
2823 theLoops.AddNewLoop();
2824 theLoops.AddEdge( myLinks[0] );
2825 if ( myLinks[0].myNode2 == myLinks[1].myNode1 )
2827 theLoops.AddEdge( myLinks[1] );
2828 theLoops.AddEdge( myLinks[2] );
2832 theLoops.AddEdge( myLinks[2] );
2833 theLoops.AddEdge( myLinks[1] );
2838 while ( !theLoops.AllEdgesUsed() )
2840 EdgeLoop& loop = theLoops.AddNewLoop();
2842 // add 1st edge to a new loop
2844 for ( i1 = theLoops.myNbLoops - 1; i1 < myLinks.size(); ++i1 )
2845 if ( theLoops.AddEdge( myLinks[i1] ))
2848 EdgePart* lastEdge = & myLinks[ i1 ];
2849 EdgePart* twinEdge = getTwin( lastEdge );
2850 const SMDS_MeshNode* firstNode = lastEdge->myNode1;
2851 const SMDS_MeshNode* lastNode = lastEdge->myNode2;
2853 do // add the rest edges
2855 theLoops.myCandidates.clear(); // edges starting at lastNode
2858 // find candidate edges
2859 for ( size_t i = i1 + 1; i < myLinks.size(); ++i )
2860 if ( myLinks[ i ].myNode1 == lastNode &&
2861 &myLinks[ i ] != twinEdge &&
2862 !theLoops.myIsUsedEdge[ i ])
2864 theLoops.myCandidates.push_back( & myLinks[ i ]);
2865 nbInternal += myLinks[ i ].IsInternal();
2868 // choose among candidates
2869 if ( theLoops.myCandidates.size() == 0 )
2871 loop.myHasPending = bool( twinEdge );
2872 lastEdge = twinEdge;
2874 else if ( theLoops.myCandidates.size() == 1 )
2876 lastEdge = theLoops.myCandidates[0];
2878 else if ( nbInternal == 1 && !lastEdge->IsInternal() )
2880 lastEdge = theLoops.myCandidates[ !theLoops.myCandidates[0]->IsInternal() ];
2884 gp_Vec lastVec = *lastEdge;
2885 double maxAngle = -2 * M_PI;
2886 for ( size_t i = 0; i < theLoops.myCandidates.size(); ++i )
2888 double angle = lastVec.AngleWithRef( *theLoops.myCandidates[i], theFaceNorm );
2889 if ( angle > maxAngle )
2892 lastEdge = theLoops.myCandidates[i];
2896 theLoops.AddEdge( *lastEdge );
2897 lastNode = lastEdge->myNode2;
2898 twinEdge = getTwin( lastEdge );
2900 while ( lastNode != firstNode );
2903 if ( twinEdge == & myLinks[ i1 ])
2904 loop.myHasPending = true;
2906 } // while ( !theLoops.AllEdgesUsed() )
2911 //================================================================================
2913 * \brief Erase loops that are cut off by face intersections
2915 //================================================================================
2917 void CutFace::CutOffLoops( EdgeLoopSet& theLoops,
2918 const double theSign,
2919 const std::vector< gp_XYZ >& theNormals,
2920 std::vector< EdgePart >& theCutOffLinks,
2921 TLinkMap& theCutOffCoplanarLinks) const
2924 boost::container::flat_set< const SMDS_MeshElement* > checkedCoplanar;
2926 for ( size_t i = 0; i < myLinks.size(); ++i )
2928 if ( !myLinks[i].myFace )
2931 EdgeLoop* loop = theLoops.GetLoopOf( & myLinks[i] );
2932 if ( !loop || loop->myLinks.empty() || loop->myHasPending )
2935 bool toErase, isCoplanar = ( myLinks[i].myIndex == EdgePart::_COPLANAR );
2937 gp_Vec iniNorm = theNormals[ myInitFace->GetID() ];
2940 toErase = ( myLinks[i].myFace->GetID() > myInitFace->GetID() );
2942 const EdgePart* twin = getTwin( & myLinks[i] );
2943 if ( !twin || twin->myFace == myLinks[i].myFace )
2945 // only one co-planar face includes myLinks[i]
2946 gp_Vec inFaceDir = iniNorm ^ myLinks[i];
2947 gp_XYZ edgePnt = SMESH_NodeXYZ( myLinks[i].myNode1 );
2948 for ( int iN = 0; iN < 3; ++iN )
2950 gp_Vec inCutFaceDir = ( SMESH_NodeXYZ( myLinks[i].myFace->GetNode( iN )) - edgePnt );
2951 if ( inCutFaceDir * inFaceDir < 0 )
2961 gp_Vec cutNorm = theNormals[ myLinks[i].myFace->GetID() ];
2962 gp_Vec inFaceDir = iniNorm ^ myLinks[i];
2964 toErase = inFaceDir * cutNorm * theSign < 0;
2967 // erase a neighboring loop
2969 if ( const EdgePart* twin = getTwin( & myLinks[i] ))
2970 loop = theLoops.GetLoopOf( twin );
2971 toErase = ( loop && !loop->myLinks.empty() );
2974 if ( toErase ) // do not erase if cutFace is connected to a co-planar cutFace
2976 checkedCoplanar.clear();
2977 for ( size_t iE = 0; iE < myLinks.size() && toErase; ++iE )
2979 if ( !myLinks[iE].myFace || myLinks[iE].myIndex != EdgePart::_COPLANAR )
2981 bool isAdded = checkedCoplanar.insert( myLinks[iE].myFace ).second;
2984 toErase = ( SMESH_MeshAlgos::NbCommonNodes( myLinks[i ].myFace,
2985 myLinks[iE].myFace ) < 1 );
2994 // remember whole sides of myInitFace that are cut off
2995 for ( size_t iE = 0; iE < loop->myLinks.size(); ++iE )
2997 if ( !loop->myLinks[ iE ]->myFace &&
2998 !loop->myLinks[ iE ]->IsInternal() )// &&
2999 // !loop->myLinks[ iE ]->myNode1->isMarked() && // cut nodes are marked
3000 // !loop->myLinks[ iE ]->myNode2->isMarked() )
3002 int i = loop->myLinks[ iE ]->myIndex;
3003 sideEdge.Set( myInitFace->GetNode ( i ),
3004 myInitFace->GetNodeWrap( i+1 ));
3005 theCutOffLinks.push_back( sideEdge );
3007 if ( !sideEdge.IsSame( *loop->myLinks[ iE ] )) // nodes replaced
3009 theCutOffLinks.push_back( *loop->myLinks[ iE ] );
3012 else if ( IsCoplanar( loop->myLinks[ iE ]))
3014 // propagate erasure to a co-planar face
3015 theCutOffLinks.push_back( *loop->myLinks[ iE ]);
3017 else if ( loop->myLinks[ iE ]->myFace &&
3018 loop->myLinks[ iE ]->IsInternal() )
3019 theCutOffLinks.push_back( *loop->myLinks[ iE ]);
3023 theLoops.Erase( loop );
3030 //================================================================================
3032 * \brief Check if the face has cut edges
3034 //================================================================================
3036 bool CutFace::IsCut() const
3038 if ( myLinks.size() > 3 )
3041 if ( myLinks.size() == 3 )
3042 for ( size_t i = 0; i < 3; ++i )
3043 if ( myLinks[i].myFace )
3049 //================================================================================
3051 * \brief Check if an edge is produced by a co-planar cut
3053 //================================================================================
3055 bool CutFace::IsCoplanar( const EdgePart* edge ) const
3057 if ( edge->myIndex == EdgePart::_COPLANAR )
3059 const EdgePart* twin = getTwin( edge );
3060 return ( !twin || twin->myIndex == EdgePart::_COPLANAR );
3065 //================================================================================
3067 * \brief Replace _COPLANAR cut edge by _INTERNAL or vice versa
3069 //================================================================================
3071 bool EdgePart::ReplaceCoplanar( const EdgePart& e )
3073 if ( myIndex + e.myIndex == _COPLANAR + _INTERNAL )
3075 //check if the faces are connected
3076 int nbCommonNodes = 0;
3077 if ( e.myFace && myFace )
3078 nbCommonNodes = SMESH_MeshAlgos::NbCommonNodes( e.myFace, myFace );
3079 bool toReplace = (( myIndex == _INTERNAL && nbCommonNodes > 1 ) ||
3080 ( myIndex == _COPLANAR && nbCommonNodes < 2 ));
3083 myIndex = e.myIndex;
3093 //================================================================================
3095 * \brief Create an offsetMesh of given faces
3096 * \param [in] faceIt - the input faces
3097 * \param [out] new2OldFaces - history of faces (new face -> old face ID)
3098 * \param [out] new2OldNodes - history of nodes (new node -> old node ID)
3099 * \return SMDS_Mesh* - the new offset mesh, a caller should delete
3101 //================================================================================
3103 SMDS_Mesh* SMESH_MeshAlgos::MakeOffset( SMDS_ElemIteratorPtr theFaceIt,
3104 SMDS_Mesh& theSrcMesh,
3105 const double theOffset,
3106 const bool theFixIntersections,
3107 TElemIntPairVec& theNew2OldFaces,
3108 TNodeIntPairVec& theNew2OldNodes)
3110 if ( theSrcMesh.GetMeshInfo().NbFaces( ORDER_QUADRATIC ) > 0 )
3111 throw SALOME_Exception( "Offset of quadratic mesh not supported" );
3112 if ( theSrcMesh.GetMeshInfo().NbFaces() > theSrcMesh.GetMeshInfo().NbTriangles() )
3113 throw SALOME_Exception( "Offset of non-triangular mesh not supported" );
3115 SMDS_Mesh* newMesh = new SMDS_Mesh;
3116 theNew2OldFaces.clear();
3117 theNew2OldNodes.clear();
3118 theNew2OldFaces.push_back
3119 ( std::make_pair(( const SMDS_MeshElement*) 0, 0)); // to have index == face->GetID()
3121 // copy input faces to the newMesh keeping IDs of nodes
3123 double minNodeDist = 1e100;
3125 std::vector< const SMDS_MeshNode* > nodes;
3126 while ( theFaceIt->more() )
3128 const SMDS_MeshElement* face = theFaceIt->next();
3129 if ( face->GetType() != SMDSAbs_Face ) continue;
3132 nodes.assign( face->begin_nodes(), face->end_nodes() );
3133 for ( size_t i = 0; i < nodes.size(); ++i )
3135 const SMDS_MeshNode* newNode = newMesh->FindNode( nodes[i]->GetID() );
3138 SMESH_NodeXYZ xyz( nodes[i] );
3139 newNode = newMesh->AddNodeWithID( xyz.X(), xyz.Y(), xyz.Z(), nodes[i]->GetID() );
3140 theNew2OldNodes.push_back( std::make_pair( newNode, nodes[i]->GetID() ));
3144 const SMDS_MeshElement* newFace = 0;
3145 switch ( face->GetEntityType() )
3147 case SMDSEntity_Triangle:
3148 newFace = newMesh->AddFace( nodes[0],nodes[1],nodes[2] );
3150 case SMDSEntity_Quad_Triangle:
3151 newFace = newMesh->AddFace( nodes[0],nodes[1],nodes[2],
3152 nodes[3],nodes[4],nodes[5] );
3154 case SMDSEntity_BiQuad_Triangle:
3155 newFace = newMesh->AddFace( nodes[0],nodes[1],nodes[2],
3156 nodes[3],nodes[4],nodes[5],nodes[6] );
3158 case SMDSEntity_Quadrangle:
3159 newFace = newMesh->AddFace( nodes[0],nodes[1],nodes[2],nodes[3] );
3161 case SMDSEntity_Quad_Quadrangle:
3162 newFace = newMesh->AddFace( nodes[0],nodes[1],nodes[2],nodes[3],
3163 nodes[4],nodes[5],nodes[6],nodes[7] );
3165 case SMDSEntity_BiQuad_Quadrangle:
3166 newFace = newMesh->AddFace( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4],
3167 nodes[5],nodes[6],nodes[7],nodes[8] );
3169 case SMDSEntity_Polygon:
3170 newFace = newMesh->AddPolygonalFace( nodes );
3172 case SMDSEntity_Quad_Polygon:
3173 newFace = newMesh->AddQuadPolygonalFace( nodes );
3178 theNew2OldFaces.push_back( std::make_pair( newFace, face->GetID() ));
3180 SMESH_NodeXYZ pPrev = nodes.back(), p;
3181 for ( size_t i = 0; i < nodes.size(); ++i )
3184 double dist = ( pPrev - p ).SquareModulus();
3185 if ( dist < minNodeDist && dist > std::numeric_limits<double>::min() )
3189 } // while ( faceIt->more() )
3192 // compute normals to faces
3193 std::vector< gp_XYZ > normals( theNew2OldFaces.size() );
3194 for ( size_t i = 1; i < normals.size(); ++i )
3196 if ( !SMESH_MeshAlgos::FaceNormal( theNew2OldFaces[i].first, normals[i] ))
3197 normals[i].SetCoord( 0,0,0 ); // TODO find norm by neighbors
3200 const double sign = ( theOffset < 0 ? -1 : +1 );
3201 const double tol = Min( 1e-3 * Sqrt( minNodeDist ),
3202 1e-2 * theOffset * sign );
3204 // translate new nodes by normal to input faces
3206 std::vector< const SMDS_MeshNode* > multiNormalNodes;
3207 for ( size_t i = 0; i < theNew2OldNodes.size(); ++i )
3209 const SMDS_MeshNode* newNode = theNew2OldNodes[i].first;
3211 if ( getTranslatedPosition( newNode, theOffset, tol*10., sign, normals, theSrcMesh, newXYZ ))
3212 newMesh->MoveNode( newNode, newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
3214 multiNormalNodes.push_back( newNode );
3216 // make multi-normal translation
3217 std::vector< SMESH_NodeXYZ > multiPos(10);
3218 for ( size_t i = 0; i < multiNormalNodes.size(); ++i )
3220 const SMDS_MeshNode* newNode = multiNormalNodes[i];
3221 newNode->setIsMarked( true );
3222 SMESH_NodeXYZ oldXYZ = newNode;
3224 for ( SMDS_ElemIteratorPtr fIt = newNode->GetInverseElementIterator(); fIt->more(); )
3226 const SMDS_MeshElement* newFace = fIt->next();
3227 const int faceIndex = newFace->GetID();
3228 const gp_XYZ& oldNorm = normals[ faceIndex ];
3229 const gp_XYZ newXYZ = oldXYZ + oldNorm * theOffset;
3230 if ( multiPos.empty() )
3232 newMesh->MoveNode( newNode, newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
3233 multiPos.emplace_back( newNode );
3238 for ( size_t iP = 0; iP < multiPos.size() && !newNode; ++iP )
3239 if (( multiPos[iP] - newXYZ ).SquareModulus() < tol * tol )
3240 newNode = multiPos[iP].Node();
3243 newNode = newMesh->AddNode( newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
3244 newNode->setIsMarked( true );
3245 theNew2OldNodes.push_back( std::make_pair( newNode, 0 ));
3246 multiPos.emplace_back( newNode );
3249 if ( newNode != oldXYZ.Node() )
3251 nodes.assign( newFace->begin_nodes(), newFace->end_nodes() );
3252 nodes[ newFace->GetNodeIndex( oldXYZ.Node() )] = newNode;
3253 newMesh->ChangeElementNodes( newFace, & nodes[0], nodes.size() );
3258 if ( !theFixIntersections )
3262 // remove new faces around concave nodes (they are marked) if the faces are inverted
3264 for ( size_t i = 0; i < theNew2OldNodes.size(); ++i )
3266 const SMDS_MeshNode* newNode = theNew2OldNodes[i].first;
3267 //const SMDS_MeshNode* oldNode = theNew2OldNodes[i].second;
3268 if ( newNode->isMarked() )
3270 //gp_XYZ moveVec = sign * ( SMESH_NodeXYZ( newNode ) - SMESH_NodeXYZ( oldNode ));
3272 //bool haveInverseFace = false;
3273 for ( SMDS_ElemIteratorPtr fIt = newNode->GetInverseElementIterator(); fIt->more(); )
3275 const SMDS_MeshElement* newFace = fIt->next();
3276 const int faceIndex = newFace->GetID();
3277 const gp_XYZ& oldNorm = normals[ faceIndex ];
3278 if ( !SMESH_MeshAlgos::FaceNormal( newFace, faceNorm, /*normalize=*/false ) ||
3279 //faceNorm * moveVec < 0 )
3280 faceNorm * oldNorm < 0 )
3282 //haveInverseFace = true;
3283 theNew2OldFaces[ faceIndex ].first = 0;
3284 newMesh->RemoveFreeElement( newFace );
3288 // if ( haveInverseFace )
3290 // newMesh->MoveNode( newNode, oldNode->X(), oldNode->Y(), oldNode->Z() );
3292 // for ( SMDS_ElemIteratorPtr fIt = newNode->GetInverseElementIterator(); fIt->more(); )
3294 // const SMDS_MeshElement* newFace = fIt->next();
3295 // if ( !SMESH_MeshAlgos::FaceNormal( newFace, normals[ newFace->GetID() ] ))
3296 // normals[i].SetCoord( 0,0,0 ); // TODO find norm by neighbors
3300 // mark all new nodes located closer than theOffset from theSrcMesh
3303 removeSmallFaces( newMesh, theNew2OldFaces, tol*tol );
3305 // ==================================================
3306 // find self-intersections of new faces and fix them
3307 // ==================================================
3309 std::unique_ptr< SMESH_ElementSearcher > fSearcher
3310 ( SMESH_MeshAlgos::GetElementSearcher( *newMesh, tol ));
3312 Intersector intersector( newMesh, tol, normals );
3314 std::vector< const SMDS_MeshElement* > closeFaces;
3315 std::vector< SMESH_NodeXYZ > faceNodes;
3318 for ( size_t iF = 1; iF < theNew2OldFaces.size(); ++iF )
3320 const SMDS_MeshElement* newFace = theNew2OldFaces[iF].first;
3321 if ( !newFace ) continue;
3322 faceNodes.assign( newFace->begin_nodes(), newFace->end_nodes() );
3324 bool isConcaveNode1 = false;
3325 for ( size_t iN = 0; iN < faceNodes.size() && !isConcaveNode1; ++iN )
3326 isConcaveNode1 = faceNodes[iN]->isMarked();
3328 // get faces close to a newFace
3331 for ( size_t i = 0; i < faceNodes.size(); ++i )
3332 faceBox.Add( faceNodes[i] );
3333 faceBox.Enlarge( tol );
3335 fSearcher->GetElementsInBox( faceBox, SMDSAbs_Face, closeFaces );
3337 // intersect the newFace with closeFaces
3339 for ( size_t i = 0; i < closeFaces.size(); ++i )
3341 const SMDS_MeshElement* closeFace = closeFaces[i];
3342 if ( closeFace->GetID() <= newFace->GetID() )
3343 continue; // this pair already treated
3345 // do not intersect connected faces if they have no concave nodes
3346 int nbCommonNodes = 0;
3347 for ( size_t iN = 0; iN < faceNodes.size(); ++iN )
3348 nbCommonNodes += ( closeFace->GetNodeIndex( faceNodes[iN].Node() ) >= 0 );
3350 if ( !isConcaveNode1 )
3352 bool isConcaveNode2 = false;
3353 for ( SMDS_ElemIteratorPtr nIt = closeFace->nodesIterator(); nIt->more(); )
3354 if (( isConcaveNode2 = nIt->next()->isMarked() ))
3357 if ( !isConcaveNode2 && nbCommonNodes > 0 )
3359 if ( normals[ newFace->GetID() ] * normals[ closeFace->GetID() ] < 1.0 )
3360 continue; // not co-planar
3364 intersector.Cut( newFace, closeFace, nbCommonNodes );
3367 intersector.MakeNewFaces( theNew2OldFaces, theNew2OldNodes, sign, /*optimize=*/true );