Salome HOME
03e695a1dee06a90d35c07f9975755e9fcc0c3d7
[modules/smesh.git] / src / SMESHUtils / SMESH_MeshAlgos.hxx
1 // Copyright (C) 2007-2019  CEA/DEN, EDF R&D, 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    * Return a face having linked nodes n1 and n2 and which is
173    * - not in avoidSet,
174    * - in elemSet provided that !elemSet.empty()
175    * i1 and i2 optionally returns indices of n1 and n2
176    */
177   SMESHUtils_EXPORT
178   const SMDS_MeshElement* FindFaceInSet(const SMDS_MeshNode*    n1,
179                                         const SMDS_MeshNode*    n2,
180                                         const TIDSortedElemSet& elemSet,
181                                         const TIDSortedElemSet& avoidSet,
182                                         int*                    i1=0,
183                                         int*                    i2=0);
184   /*!
185    * \brief Calculate normal of a mesh face
186    */
187   SMESHUtils_EXPORT
188   bool FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool normalized=true);
189
190   /*!
191    * \brief Return nodes common to two elements
192    */
193   SMESHUtils_EXPORT
194   std::vector< const SMDS_MeshNode*> GetCommonNodes(const SMDS_MeshElement* e1,
195                                                     const SMDS_MeshElement* e2);
196   /*!
197    * \brief Return true if node1 encounters first in the face and node2, after.
198    *        The nodes are supposed to be neighbor nodes in the face.
199    */
200   SMESHUtils_EXPORT
201   bool IsRightOrder( const SMDS_MeshElement* face,
202                      const SMDS_MeshNode*    node0,
203                      const SMDS_MeshNode*    node1 );
204
205   typedef std::vector< std::vector< const SMDS_MeshElement* > > TElemGroupVector;
206   typedef std::vector< std::vector< const SMDS_MeshNode* > >    TNodeGroupVector;
207   /*!
208    * \brief Partition given 1D elements into groups of contiguous edges.
209    *        A node where number of meeting edges != 2 is a group end.
210    *        An optional startNode is used to orient groups it belongs to.
211    * \return a list of edge groups and a list of corresponding node groups.
212    *         If a group is closed, the first and last nodes of the group are same.
213    */
214   SMESHUtils_EXPORT
215   void Get1DBranches( SMDS_ElemIteratorPtr edgeIt,
216                       TElemGroupVector&    edgeGroups,
217                       TNodeGroupVector&    nodeGroups,
218                       const SMDS_MeshNode* startNode = 0 );
219
220   /*!
221    * \brief Mark elements given by SMDS_Iterator
222    */
223   template< class ElemIter >
224   void MarkElems( ElemIter it, const bool isMarked )
225   {
226     while ( it->more() ) it->next()->setIsMarked( isMarked );
227   }
228   /*!
229    * \brief Mark elements given by std iterators
230    */
231   template< class ElemIter >
232   void MarkElems( ElemIter it, ElemIter end, const bool isMarked )
233   {
234     for ( ; it != end; ++it ) (*it)->setIsMarked( isMarked );
235   }
236   /*!
237    * \brief Mark nodes of elements given by SMDS_Iterator
238    */
239   template< class ElemIter >
240   void MarkElemNodes( ElemIter it, const bool isMarked, const bool markElem = false )
241   {
242     if ( markElem )
243       while ( it->more() ) {
244         const SMDS_MeshElement* e = it->next();
245         e->setIsMarked( isMarked );
246         MarkElems( e->nodesIterator(), isMarked );
247       }
248     else
249       while ( it->more() )
250         MarkElems( it->next()->nodesIterator(), isMarked );
251   }
252   /*!
253    * \brief Mark elements given by std iterators
254    */
255   template< class ElemIter >
256   void MarkElemNodes( ElemIter it, ElemIter end, const bool isMarked, const bool markElem = false )
257   {
258     if ( markElem )
259       for ( ; it != end; ++it ) {
260         (*it)->setIsMarked( isMarked );
261         MarkElems( (*it)->nodesIterator(), isMarked );
262       }
263     else
264       for ( ; it != end; ++it )
265         MarkElems( (*it)->nodesIterator(), isMarked );
266   }
267
268   // 2 nodes + optional medium node
269   struct Edge
270   {
271     const SMDS_MeshNode* _node1;
272     const SMDS_MeshNode* _node2;
273     const SMDS_MeshNode* _medium;
274   };
275
276   /*!
277    * Return sharp edges of faces and non-manifold ones.
278    * Optionally adds existing edges to the result. Angle is in degrees.
279    */
280   SMESHUtils_EXPORT
281   std::vector< Edge > FindSharpEdges( SMDS_Mesh* mesh,
282                                       double     angle,
283                                       bool       addExisting );
284
285   /*!
286    * Distribute all faces of the mesh between groups using given edges.
287    */
288   SMESHUtils_EXPORT
289   std::vector< std::vector< const SMDS_MeshElement* > >
290   SeparateFacesByEdges( SMDS_Mesh* mesh, const std::vector< Edge >& edges );
291
292
293   typedef std::vector<const SMDS_MeshNode*> TFreeBorder;
294   typedef std::vector<TFreeBorder>          TFreeBorderVec;
295   struct TFreeBorderPart
296   {
297     int _border; // border index within a TFreeBorderVec
298     int _node1;  // node index within the border-th TFreeBorder
299     int _node2;
300     int _nodeLast;
301   };
302   typedef std::vector<TFreeBorderPart>  TCoincidentGroup;
303   typedef std::vector<TCoincidentGroup> TCoincidentGroupVec;
304   struct CoincidentFreeBorders
305   {
306     TFreeBorderVec      _borders;          // nodes of all free borders
307     TCoincidentGroupVec _coincidentGroups; // groups of coincident parts of borders
308   };
309
310   /*!
311    * Returns TFreeBorder's coincident within the given tolerance.
312    * If the tolerance <= 0.0 then one tenth of an average size of elements adjacent
313    * to free borders being compared is used.
314    */
315   SMESHUtils_EXPORT
316   void FindCoincidentFreeBorders(SMDS_Mesh&              mesh,
317                                  double                  tolerance,
318                                  CoincidentFreeBorders & foundFreeBordes);
319   // Implemented in ./SMESH_FreeBorders.cxx
320
321   /*!
322    * Returns all or only closed TFreeBorder's.
323    * Optionally check if the mesh is manifold and if faces are correctly oriented.
324    */
325   SMESHUtils_EXPORT
326   void FindFreeBorders(SMDS_Mesh&       mesh,
327                        TFreeBorderVec & foundFreeBordes,
328                        const bool       closedOnly,
329                        bool*            isManifold = 0,
330                        bool*            isGoodOri = 0);
331   // Implemented in ./SMESH_FreeBorders.cxx
332
333   /*!
334    * Fill a hole defined by a TFreeBorder with 2D elements.
335    */
336   SMESHUtils_EXPORT
337   void FillHole(const TFreeBorder &                   freeBorder,
338                 SMDS_Mesh&                            mesh,
339                 std::vector<const SMDS_MeshElement*>& newFaces);
340   // Implemented in ./SMESH_FillHole.cxx
341
342   /*!
343    * \brief Find nodes whose merge makes the element invalid
344    */
345   SMESHUtils_EXPORT
346   void DeMerge(const SMDS_MeshElement*              elem,
347                std::vector< const SMDS_MeshNode* >& newNodes,
348                std::vector< const SMDS_MeshNode* >& noMergeNodes);
349   // Implemented in SMESH_DeMerge.cxx
350
351
352   typedef std::vector< std::pair< const SMDS_MeshElement*, int > > TElemIntPairVec;
353   typedef std::vector< std::pair< const SMDS_MeshNode*,    int > > TNodeIntPairVec;
354   /*!
355    * \brief Create an offset mesh of given faces
356    *  \param [in] faceIt - the input faces
357    *  \param [in] theFixIntersections - to fix self intersections of the offset mesh or not
358    *  \param [out] new2OldFaces - history of faces
359    *  \param [out] new2OldNodes - history of nodes
360    *  \return SMDS_Mesh* - the new offset mesh, a caller should delete
361    */
362   SMESHUtils_EXPORT
363   SMDS_Mesh* MakeOffset( SMDS_ElemIteratorPtr faceIt,
364                          SMDS_Mesh&           mesh,
365                          const double         offset,
366                          const bool           theFixIntersections,
367                          TElemIntPairVec&           new2OldFaces,
368                          TNodeIntPairVec&           new2OldNodes );
369   // Implemented in ./SMESH_Offset.cxx
370
371
372   //=======================================================================
373   /*!
374    * \brief Cut faces of a triangular mesh.
375    *  Usage work-flow: 1) call Cut() methods as many times as needed
376    *                   2) call MakeNewFaces() to really modify the mesh faces
377    */
378   //=======================================================================
379   // implemented in SMESH_Offset.cxx
380
381   class SMESHUtils_EXPORT Intersector
382   {
383   public:
384     Intersector( SMDS_Mesh* mesh, double tol, const std::vector< gp_XYZ >& normals );
385     ~Intersector();
386
387     //! Compute cut of two faces of the mesh
388     void Cut( const SMDS_MeshElement* face1,
389               const SMDS_MeshElement* face2,
390               const int               nbCommonNodes = -1 );
391
392     //! Store a face cut by a line given by its ends lying either on face edges or inside the face.
393     //  Line ends are accompanied by indices of intersected face edges.
394     //  Edge index is <0 if a line end is inside the face.
395     void Cut( const SMDS_MeshElement* face,
396               SMESH_NodeXYZ&          lineEnd1,
397               int                     edgeIndex1,
398               SMESH_NodeXYZ&          lineEnd2,
399               int                     edgeIndex2 );
400
401     //! Split all faces intersected by Cut() methods.
402     //  theSign = (-1|1) is used to choose which part of a face cut by another one to remove.
403     //  1 means to remove a part opposite to face normal.
404     //  Optionally optimize quality of split faces by edge swapping.
405     void MakeNewFaces( SMESH_MeshAlgos::TElemIntPairVec& theNew2OldFaces,
406                        SMESH_MeshAlgos::TNodeIntPairVec& theNew2OldNodes,
407                        const double                      theSign = 1.,
408                        const bool                        theOptimize = false );
409
410     typedef std::vector< SMESH_NodeXYZ > TFace;
411
412     //! Cut a face by planes, whose normals point to parts to keep.
413     //  Return true if the whole face is cut off
414     static bool CutByPlanes(const SMDS_MeshElement*       face,
415                             const std::vector< gp_Ax1 > & planes,
416                             const double                  tol,
417                             std::vector< TFace > &        newFaceConnectivity );
418
419   private:
420     struct Algo;
421     Algo* myAlgo;
422   };
423
424   //=======================================================================
425   /*!
426    * \brief Divide a mesh face into triangles
427    */
428   //=======================================================================
429   // Implemented in ./SMESH_Triangulate.cxx
430
431   class SMESHUtils_EXPORT Triangulate
432   {
433   public:
434
435     Triangulate(bool optimize=false);
436     ~Triangulate();
437
438     static int GetNbTriangles( const SMDS_MeshElement* face );
439
440     int GetTriangles( const SMDS_MeshElement*             face,
441                       std::vector< const SMDS_MeshNode*>& nodes);
442   private:
443
444     bool triangulate( std::vector< const SMDS_MeshNode*>& nodes, const size_t nbNodes );
445
446     /*!
447      * \brief Vertex of a polygon. Together with 2 neighbor Vertices represents a triangle
448      */
449     struct PolyVertex
450     {
451       SMESH_NodeXYZ _nxyz;
452       size_t        _index;
453       gp_XY         _xy;
454       PolyVertex*   _prev;
455       PolyVertex*   _next;
456
457       void   SetNodeAndNext( const SMDS_MeshNode* n, PolyVertex& v, size_t index );
458       void   GetTriaNodes( const SMDS_MeshNode** nodes, size_t* nodeIndices) const;
459       double TriaArea() const;
460       bool   IsInsideTria( const PolyVertex* v );
461       PolyVertex* Delete();
462     };
463     struct Optimizer;
464
465     std::vector< PolyVertex > _pv;
466     std::vector< size_t >     _nodeIndex;
467     Optimizer*                _optimizer;
468   };
469
470   // structure used in MakePolyLine() to define a cutting plane
471   struct PolySegment
472   {
473     // 2 points, each defined as follows:
474     // ( myNode1 &&  myNode2 ) ==> point is in the middle of an edge defined by two nodes
475     // ( myNode1 && !myNode2 ) ==> point is at myNode1 of a some face
476     // else                    ==> point is at myXYZ
477     const SMDS_MeshNode*    myNode1[2];
478     const SMDS_MeshNode*    myNode2[2];
479     gp_XYZ                  myXYZ  [2];
480
481     // face on which myXYZ projects (found by MakePolyLine())
482     const SMDS_MeshElement* myFace [2];
483
484     // vector on the plane; to use a default plane set vector = (0,0,0)
485     gp_Vec myVector;
486
487     // point returning coordinates of a middle of the two points, projected to mesh
488     gp_Pnt myMidProjPoint;
489   };
490   typedef std::vector<PolySegment> TListOfPolySegments;
491
492   /*!
493    * \brief Create a polyline consisting of 1D mesh elements each lying on a 2D element of
494    *        the initial mesh. Positions of new nodes are found by cutting the mesh by the
495    *        plane passing through pairs of points specified by each PolySegment structure.
496    *        If there are several paths connecting a pair of points, the shortest path is
497    *        selected by the module. Position of the cutting plane is defined by the two
498    *        points and an optional vector lying on the plane specified by a PolySegment.
499    *        By default the vector is defined by Mesh module as following. A middle point
500    *        of the two given points is computed. The middle point is projected to the mesh.
501    *        The vector goes from the middle point to the projection point. In case of planar
502    *        mesh, the vector is normal to the mesh.
503    *  \param [inout] segments - PolySegment's defining positions of cutting planes.
504    *        Return the used vector and position of the middle point.
505    *  \param [in] group - an optional group where created mesh segments will
506    *        be added.
507    */
508   // Implemented in ./SMESH_PolyLine.cxx
509   SMESHUtils_EXPORT
510   void MakePolyLine( SMDS_Mesh*                            mesh,
511                      TListOfPolySegments&                  segments,
512                      std::vector<const SMDS_MeshElement*>& newEdges,
513                      std::vector<const SMDS_MeshNode*>&    newNodes,
514                      SMDS_MeshGroup*                       group=0,
515                      SMESH_ElementSearcher*                searcher=0);
516
517   /*!
518    * Create a slot of given width around given 1D elements lying on a triangle mesh.
519    * The slot is consrtucted by cutting faces by cylindrical surfaces made around each segment.
520    * \return Edges located at the slot boundary
521    */
522   // Implemented in ./SMESH_Slot.cxx
523   SMESHUtils_EXPORT
524   std::vector< Edge > MakeSlot( SMDS_ElemIteratorPtr segmentIt,
525                                 double               width,
526                                 SMDS_Mesh*           mesh);
527
528 } // namespace SMESH_MeshAlgos
529
530 #endif