Salome HOME
Update of CheckDone
[modules/smesh.git] / src / SMESHUtils / SMESH_MeshAlgos.hxx
1 // Copyright (C) 2007-2024  CEA, EDF, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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.
10 //
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.
15 //
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
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 // File      : SMESH_MeshAlgos.hxx
23 // Created   : Tue Apr 30 18:00:36 2013
24 // Author    : Edward AGAPOV (eap)
25
26 // Initially this file held some low level algorithms extracted from SMESH_MeshEditor
27 // to make them accessible from Controls package, and more
28
29
30 #ifndef __SMESH_MeshAlgos_HXX__
31 #define __SMESH_MeshAlgos_HXX__
32
33 #include "SMESH_Utils.hxx"
34
35 #include "SMDSAbs_ElementType.hxx"
36 #include "SMDS_ElemIterator.hxx"
37 #include "SMESH_TypeDefs.hxx"
38
39 #include <TopAbs_State.hxx>
40 #include <gp_Pnt.hxx>
41 #include <gp_Vec.hxx>
42
43 #include <vector>
44
45 class Bnd_B3d;
46 class gp_Ax1;
47 class SMDS_Mesh;
48 class SMDS_MeshElement;
49 class SMDS_MeshGroup;
50 class SMDS_MeshNode;
51
52 //=======================================================================
53 /*!
54  * \brief Searcher for the node closest to a point
55  */
56 //=======================================================================
57
58 struct SMESHUtils_EXPORT SMESH_NodeSearcher
59 {
60   virtual const SMDS_MeshNode* FindClosestTo( const gp_Pnt& pnt ) = 0;
61   virtual void MoveNode( const SMDS_MeshNode* node, const gp_Pnt& toPnt ) = 0;
62   virtual int  FindNearPoint(const gp_Pnt&                        point,
63                              const double                         tolerance,
64                              std::vector< const SMDS_MeshNode* >& foundNodes) = 0;
65   virtual ~SMESH_NodeSearcher() {}
66 };
67
68 //=======================================================================
69 /*!
70  * \brief Searcher for elements
71  */
72 //=======================================================================
73
74 struct SMESHUtils_EXPORT SMESH_ElementSearcher
75 {
76   /*!
77    * \brief Find elements of given type where the given point is IN or ON.
78    *        Returns nb of found elements and elements them-selves.
79    *
80    * 'ALL' type means elements of any type excluding nodes and 0D elements
81    */
82   virtual int FindElementsByPoint(const gp_Pnt&                           point,
83                                   SMDSAbs_ElementType                     type,
84                                   std::vector< const SMDS_MeshElement* >& foundElems) = 0;
85   /*!
86    * \brief Return an element most close to the given point
87    */
88   virtual const SMDS_MeshElement* FindClosestTo( const gp_Pnt&       point,
89                                                  SMDSAbs_ElementType type) = 0;
90   /*!
91    * \brief Return elements possibly intersecting the line
92    */
93   virtual void GetElementsNearLine( const gp_Ax1&                           line,
94                                     SMDSAbs_ElementType                     type,
95                                     std::vector< const SMDS_MeshElement* >& foundElems) = 0;
96   /*!
97    * \brief Return elements whose bounding box intersects a sphere
98    */
99   virtual void GetElementsInSphere( const gp_XYZ&                           center,
100                                     const double                            radius,
101                                     SMDSAbs_ElementType                     type,
102                                     std::vector< const SMDS_MeshElement* >& foundElems) = 0;
103   /*!
104    * \brief Return elements whose bounding box intersects a given bounding box
105    */
106   virtual void GetElementsInBox( const Bnd_B3d&                          box,
107                                  SMDSAbs_ElementType                     type,
108                                  std::vector< const SMDS_MeshElement* >& foundElems) = 0;
109   /*!
110    * \brief Find out if the given point is out of closed 2D mesh.
111    */
112   virtual TopAbs_State GetPointState(const gp_Pnt& point) = 0;
113
114   /*!
115    * \brief Return a projection of a given point to a 2D mesh.
116    *        Optionally return the closest face
117    */
118   virtual gp_XYZ Project(const gp_Pnt&            point,
119                          SMDSAbs_ElementType      type,
120                          const SMDS_MeshElement** closestFace= 0) = 0;
121
122   virtual ~SMESH_ElementSearcher();
123 };
124
125 namespace SMESH_MeshAlgos
126 {
127   /*!
128    * \brief Return SMESH_NodeSearcher. The caller is responsible for deleting it
129    */
130   SMESHUtils_EXPORT
131   SMESH_NodeSearcher* GetNodeSearcher( SMDS_Mesh& mesh );
132
133   SMESHUtils_EXPORT
134   SMESH_NodeSearcher* GetNodeSearcher( SMDS_ElemIteratorPtr elemIt );
135
136   /*!
137    * \brief Return SMESH_ElementSearcher. The caller is responsible for deleting it
138    */
139   SMESHUtils_EXPORT
140   SMESH_ElementSearcher* GetElementSearcher( SMDS_Mesh& mesh,
141                                              double     tolerance=-1.);
142   SMESHUtils_EXPORT
143   SMESH_ElementSearcher* GetElementSearcher( SMDS_Mesh& mesh,
144                                              SMDS_ElemIteratorPtr elemIt,
145                                              double     tolerance=-1. );
146
147
148   /*!
149    * \brief Return true if the point is IN or ON of the element
150    */
151   SMESHUtils_EXPORT
152   bool IsOut( const SMDS_MeshElement* element, const gp_Pnt& point, double tol );
153
154   SMESHUtils_EXPORT
155   double GetDistance( const SMDS_MeshElement* elem, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
156
157   SMESHUtils_EXPORT
158   double GetDistance( const SMDS_MeshEdge* edge, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
159
160   SMESHUtils_EXPORT
161   double GetDistance( const SMDS_MeshFace* face, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
162
163   SMESHUtils_EXPORT
164   double GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
165
166   SMESHUtils_EXPORT
167   void GetBarycentricCoords( const gp_XY& point,
168                              const gp_XY& t0, const gp_XY& t1, const gp_XY& t2,
169                              double &    bc0, double &    bc1);
170
171   /*!
172    * \brief Intersect volume by a ray
173    */
174   bool IntersectRayVolume( const gp_Ax1& ray, const double rayLen,
175                            const SMDS_MeshElement* vol,
176                            double & tMin, double & tMax,
177                            int & iFacetMin, int & iFacetMax);
178   /*!
179    * Return a face having linked nodes n1 and n2 and which is
180    * - not in avoidSet,
181    * - in elemSet provided that !elemSet.empty()
182    * i1 and i2 optionally returns indices of n1 and n2
183    */
184   SMESHUtils_EXPORT
185   const SMDS_MeshElement* FindFaceInSet(const SMDS_MeshNode*    n1,
186                                         const SMDS_MeshNode*    n2,
187                                         const TIDSortedElemSet& elemSet,
188                                         const TIDSortedElemSet& avoidSet,
189                                         int*                    i1=0,
190                                         int*                    i2=0);
191   /*!
192    * \brief Calculate normal of a mesh face
193    */
194   SMESHUtils_EXPORT
195   bool FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool normalized=true);
196
197   /*!
198    * \brief Return number of nodes common to two elements
199    */
200   SMESHUtils_EXPORT
201   int NbCommonNodes(const SMDS_MeshElement* e1,
202                     const SMDS_MeshElement* e2);
203   /*!
204    * \brief Return nodes common to two elements
205    */
206   SMESHUtils_EXPORT
207   std::vector< const SMDS_MeshNode*> GetCommonNodes(const SMDS_MeshElement* e1,
208                                                     const SMDS_MeshElement* e2);
209   /*!
210    * \brief Return true if a node is on a boundary of 2D mesh.
211    *        Optionally returns two neighboring boundary nodes (or more in non-manifold mesh)
212    */
213   SMESHUtils_EXPORT bool IsOn2DBoundary( const SMDS_MeshNode* node,
214                                          std::vector< const SMDS_MeshNode*> * neibors = nullptr );
215   /*!
216    * \brief Return true if node1 encounters first in the face and node2, after.
217    *        The nodes are supposed to be neighbor nodes in the face.
218    */
219   SMESHUtils_EXPORT
220   bool IsRightOrder( const SMDS_MeshElement* face,
221                      const SMDS_MeshNode*    node0,
222                      const SMDS_MeshNode*    node1 );
223
224   typedef std::vector< std::vector< const SMDS_MeshElement* > > TElemGroupVector;
225   typedef std::vector< std::vector< const SMDS_MeshNode* > >    TNodeGroupVector;
226   /*!
227    * \brief Partition given 1D elements into groups of contiguous edges.
228    *        A node where number of meeting edges != 2 is a group end.
229    *        An optional startNode is used to orient groups it belongs to.
230    * \return a list of edge groups and a list of corresponding node groups.
231    *         If a group is closed, the first and last nodes of the group are same.
232    */
233   SMESHUtils_EXPORT
234   void Get1DBranches( SMDS_ElemIteratorPtr edgeIt,
235                       TElemGroupVector&    edgeGroups,
236                       TNodeGroupVector&    nodeGroups,
237                       const SMDS_MeshNode* startNode = 0 );
238
239   /*!
240    * \brief Mark elements given by SMDS_Iterator
241    * \sa SMDS_Mesh::SetAllNodesNotMarked() and SMDS_Mesh::SetAllCellsNotMarked()
242    */
243   template< class ElemIter >
244   void MarkElems( ElemIter it, const bool isMarked )
245   {
246     while ( it->more() ) it->next()->setIsMarked( isMarked );
247   }
248   /*!
249    * \brief Mark elements given by std iterators
250    */
251   template< class ElemIter >
252   void MarkElems( ElemIter it, ElemIter end, const bool isMarked )
253   {
254     for ( ; it != end; ++it ) (*it)->setIsMarked( isMarked );
255   }
256   /*!
257    * \brief Mark nodes of elements given by SMDS_Iterator
258    */
259   template< class ElemIter >
260   void MarkElemNodes( ElemIter it, const bool isMarked, const bool markElem = false )
261   {
262     if ( markElem )
263       while ( it->more() ) {
264         const SMDS_MeshElement* e = it->next();
265         e->setIsMarked( isMarked );
266         MarkElems( e->nodesIterator(), isMarked );
267       }
268     else
269       while ( it->more() )
270         MarkElems( it->next()->nodesIterator(), isMarked );
271   }
272   /*!
273    * \brief Mark elements given by std iterators
274    */
275   template< class ElemIter >
276   void MarkElemNodes( ElemIter it, ElemIter end, const bool isMarked, const bool markElem = false )
277   {
278     if ( markElem )
279       for ( ; it != end; ++it ) {
280         (*it)->setIsMarked( isMarked );
281         MarkElems( (*it)->nodesIterator(), isMarked );
282       }
283     else
284       for ( ; it != end; ++it )
285         MarkElems( (*it)->nodesIterator(), isMarked );
286   }
287
288   // 2 nodes + optional medium node
289   struct Edge
290   {
291     const SMDS_MeshNode* _node1;
292     const SMDS_MeshNode* _node2;
293     const SMDS_MeshNode* _medium;
294   };
295
296   /*!
297    * Return sharp edges of faces and non-manifold ones.
298    * Optionally adds existing edges to the result. Angle is in degrees.
299    */
300   SMESHUtils_EXPORT
301   std::vector< Edge > FindSharpEdges( SMDS_Mesh* mesh,
302                                       double     angle,
303                                       bool       addExisting );
304
305   /*!
306    * Distribute all faces of the mesh between groups using given edges.
307    */
308   SMESHUtils_EXPORT
309   std::vector< std::vector< const SMDS_MeshElement* > >
310   SeparateFacesByEdges( SMDS_Mesh* mesh, const std::vector< Edge >& edges );
311
312
313   typedef std::vector<const SMDS_MeshNode*> TFreeBorder;
314   typedef std::vector<TFreeBorder>          TFreeBorderVec;
315   struct TFreeBorderPart
316   {
317     int _border; // border index within a TFreeBorderVec
318     int _node1;  // node index within the border-th TFreeBorder
319     int _node2;
320     int _nodeLast;
321   };
322   typedef std::vector<TFreeBorderPart>  TCoincidentGroup;
323   typedef std::vector<TCoincidentGroup> TCoincidentGroupVec;
324   struct CoincidentFreeBorders
325   {
326     TFreeBorderVec      _borders;          // nodes of all free borders
327     TCoincidentGroupVec _coincidentGroups; // groups of coincident parts of borders
328   };
329
330   /*!
331    * Returns TFreeBorder's coincident within the given tolerance.
332    * If the tolerance <= 0.0 then one tenth of an average size of elements adjacent
333    * to free borders being compared is used.
334    */
335   SMESHUtils_EXPORT
336   void FindCoincidentFreeBorders(SMDS_Mesh&              mesh,
337                                  double                  tolerance,
338                                  CoincidentFreeBorders & foundFreeBordes);
339   // Implemented in ./SMESH_FreeBorders.cxx
340
341   /*!
342    * Returns all or only closed TFreeBorder's.
343    * Optionally check if the mesh is manifold and if faces are correctly oriented.
344    */
345   SMESHUtils_EXPORT
346   void FindFreeBorders(SMDS_Mesh&       mesh,
347                        TFreeBorderVec & foundFreeBordes,
348                        const bool       closedOnly,
349                        bool*            isManifold = 0,
350                        bool*            isGoodOri = 0);
351   // Implemented in ./SMESH_FreeBorders.cxx
352
353   /*!
354    * Fill a hole defined by a TFreeBorder with 2D elements.
355    */
356   SMESHUtils_EXPORT
357   void FillHole(const TFreeBorder &                   freeBorder,
358                 SMDS_Mesh&                            mesh,
359                 std::vector<const SMDS_MeshElement*>& newFaces);
360   // Implemented in ./SMESH_FillHole.cxx
361
362   /*!
363    * \brief Find nodes whose merge makes the element invalid
364    */
365   SMESHUtils_EXPORT
366   void DeMerge(const SMDS_MeshElement*              elem,
367                std::vector< const SMDS_MeshNode* >& newNodes,
368                std::vector< const SMDS_MeshNode* >& noMergeNodes);
369   // Implemented in SMESH_DeMerge.cxx
370
371
372   typedef std::vector< std::pair< const SMDS_MeshElement*, int > > TElemIntPairVec;
373   typedef std::vector< std::pair< const SMDS_MeshNode*,    int > > TNodeIntPairVec;
374   /*!
375    * \brief Create an offset mesh of given faces
376    *  \param [in] faceIt - the input faces
377    *  \param [in] theFixIntersections - to fix self intersections of the offset mesh or not
378    *  \param [out] new2OldFaces - history of faces
379    *  \param [out] new2OldNodes - history of nodes
380    *  \return SMDS_Mesh* - the new offset mesh, a caller should delete
381    */
382   SMESHUtils_EXPORT
383   SMDS_Mesh* MakeOffset( SMDS_ElemIteratorPtr faceIt,
384                          SMDS_Mesh&           mesh,
385                          const double         offset,
386                          const bool           theFixIntersections,
387                          TElemIntPairVec&           new2OldFaces,
388                          TNodeIntPairVec&           new2OldNodes );
389   // Implemented in ./SMESH_Offset.cxx
390
391
392   //=======================================================================
393   /*!
394    * \brief Cut faces of a triangular mesh.
395    *  Usage work-flow: 1) call Cut() methods as many times as needed
396    *                   2) call MakeNewFaces() to really modify the mesh faces
397    */
398   //=======================================================================
399   // implemented in SMESH_Offset.cxx
400
401   class SMESHUtils_EXPORT Intersector
402   {
403   public:
404     Intersector( SMDS_Mesh* mesh, double tol, const std::vector< gp_XYZ >& normals );
405     ~Intersector();
406
407     //! Compute cut of two faces of the mesh
408     void Cut( const SMDS_MeshElement* face1,
409               const SMDS_MeshElement* face2,
410               const int               nbCommonNodes = -1 );
411
412     //! Store a face cut by a line given by its ends lying either on face edges or inside the face.
413     //  Line ends are accompanied by indices of intersected face edges.
414     //  Edge index is <0 if a line end is inside the face.
415     void Cut( const SMDS_MeshElement* face,
416               SMESH_NodeXYZ&          lineEnd1,
417               int                     edgeIndex1,
418               SMESH_NodeXYZ&          lineEnd2,
419               int                     edgeIndex2 );
420
421     //! Split all faces intersected by Cut() methods.
422     //  theSign = (-1|1) is used to choose which part of a face cut by another one to remove.
423     //  1 means to remove a part opposite to face normal.
424     //  Optionally optimize quality of split faces by edge swapping.
425     void MakeNewFaces( SMESH_MeshAlgos::TElemIntPairVec& theNew2OldFaces,
426                        SMESH_MeshAlgos::TNodeIntPairVec& theNew2OldNodes,
427                        const double                      theSign = 1.,
428                        const bool                        theOptimize = false );
429
430     typedef std::vector< SMESH_NodeXYZ > TFace;
431
432     //! Cut a face by planes, whose normals point to parts to keep.
433     //  Return true if the whole face is cut off
434     static bool CutByPlanes(const SMDS_MeshElement*       face,
435                             const std::vector< gp_Ax1 > & planes,
436                             const double                  tol,
437                             std::vector< TFace > &        newFaceConnectivity );
438
439   private:
440     struct Algo;
441     Algo* myAlgo;
442   };
443
444   //=======================================================================
445   /*!
446    * \brief Divide a mesh face into triangles
447    */
448   //=======================================================================
449   // Implemented in ./SMESH_Triangulate.cxx
450
451   class SMESHUtils_EXPORT Triangulate
452   {
453   public:
454
455     Triangulate(bool optimize=false);
456     ~Triangulate();
457
458     static int GetNbTriangles( const SMDS_MeshElement* face );
459
460     int GetTriangles( const SMDS_MeshElement*             face,
461                       std::vector< const SMDS_MeshNode*>& nodes);
462   private:
463
464     bool triangulate( std::vector< const SMDS_MeshNode*>& nodes, const size_t nbNodes );
465
466     struct PolyVertex;
467     struct Optimizer;
468     struct Data;
469
470     Data*      _data;
471     Optimizer* _optimizer;
472   };
473
474   // structure used in MakePolyLine() to define a cutting plane
475   struct PolySegment
476   {
477     // 2 points, each defined as follows:
478     // ( myNode1 &&  myNode2 ) ==> point is in the middle of an edge defined by two nodes
479     // ( myNode1 && !myNode2 ) ==> point is at myNode1 of a some face
480     // else                    ==> point is at myXYZ
481     const SMDS_MeshNode*    myNode1[2];
482     const SMDS_MeshNode*    myNode2[2];
483     gp_XYZ                  myXYZ  [2];
484
485     // face on which myXYZ projects (found by MakePolyLine())
486     const SMDS_MeshElement* myFace [2];
487
488     // vector on the plane; to use a default plane set vector = (0,0,0)
489     gp_Vec myVector;
490
491     // point returning coordinates of a middle of the two points, projected to mesh
492     gp_Pnt myMidProjPoint;
493   };
494   typedef std::vector<PolySegment> TListOfPolySegments;
495
496   /*!
497    * \brief Create a polyline consisting of 1D mesh elements each lying on a 2D element of
498    *        the initial mesh. Positions of new nodes are found by cutting the mesh by the
499    *        plane passing through pairs of points specified by each PolySegment structure.
500    *        If there are several paths connecting a pair of points, the shortest path is
501    *        selected by the module. Position of the cutting plane is defined by the two
502    *        points and an optional vector lying on the plane specified by a PolySegment.
503    *        By default the vector is defined by Mesh module as following. A middle point
504    *        of the two given points is computed. The middle point is projected to the mesh.
505    *        The vector goes from the middle point to the projection point. In case of planar
506    *        mesh, the vector is normal to the mesh.
507    *  \param [inout] segments - PolySegment's defining positions of cutting planes.
508    *        Return the used vector and position of the middle point.
509    *  \param [in] group - an optional group where created mesh segments will
510    *        be added.
511    */
512   // Implemented in ./SMESH_PolyLine.cxx
513   SMESHUtils_EXPORT
514   void MakePolyLine( SMDS_Mesh*                            mesh,
515                      TListOfPolySegments&                  segments,
516                      std::vector<const SMDS_MeshElement*>& newEdges,
517                      std::vector<const SMDS_MeshNode*>&    newNodes,
518                      SMDS_MeshGroup*                       group=0,
519                      SMESH_ElementSearcher*                searcher=0);
520
521   /*!
522    * Create a slot of given width around given 1D elements lying on a triangle mesh.
523    * The slot is constructed by cutting faces by cylindrical surfaces made around each segment.
524    * \return Edges located at the slot boundary
525    */
526   // Implemented in ./SMESH_Slot.cxx
527   SMESHUtils_EXPORT
528   std::vector< Edge > MakeSlot( SMDS_ElemIteratorPtr             segmentIt,
529                                 double                           width,
530                                 SMDS_Mesh*                       mesh,
531                                 std::vector< SMDS_MeshGroup* > & groupsToUpdate);
532
533 } // namespace SMESH_MeshAlgos
534
535 #endif