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