Salome HOME
6a507f1b1d6fa67604d5b1cea4bde94a93693e0f
[modules/smesh.git] / src / SMESHUtils / SMESH_MeshAlgos.hxx
1 // Copyright (C) 2007-2016  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 // This file holds 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 <vector>
41
42 class gp_Pnt;
43 class gp_Ax1;
44 class Bnd_B3d;
45 class SMDS_MeshNode;
46 class SMDS_MeshElement;
47 class SMDS_Mesh;
48
49 //=======================================================================
50 /*!
51  * \brief Searcher for the node closest to a point
52  */
53 //=======================================================================
54
55 struct SMESHUtils_EXPORT SMESH_NodeSearcher
56 {
57   virtual const SMDS_MeshNode* FindClosestTo( const gp_Pnt& pnt ) = 0;
58   virtual void MoveNode( const SMDS_MeshNode* node, const gp_Pnt& toPnt ) = 0;
59   virtual int  FindNearPoint(const gp_Pnt&                        point,
60                              const double                         tolerance,
61                              std::vector< const SMDS_MeshNode* >& foundNodes) = 0;
62   virtual ~SMESH_NodeSearcher() {}
63 };
64
65 //=======================================================================
66 /*!
67  * \brief Searcher for elements
68  */
69 //=======================================================================
70
71 struct SMESHUtils_EXPORT SMESH_ElementSearcher
72 {
73   /*!
74    * \brief Find elements of given type where the given point is IN or ON.
75    *        Returns nb of found elements and elements them-selves.
76    *
77    * 'ALL' type means elements of any type excluding nodes and 0D elements
78    */
79   virtual int FindElementsByPoint(const gp_Pnt&                           point,
80                                   SMDSAbs_ElementType                     type,
81                                   std::vector< const SMDS_MeshElement* >& foundElems) = 0;
82   /*!
83    * \brief Return an element most close to the given point
84    */
85   virtual const SMDS_MeshElement* FindClosestTo( const gp_Pnt&       point,
86                                                  SMDSAbs_ElementType type) = 0;
87   /*!
88    * \brief Return elements possibly intersecting the line
89    */
90   virtual void GetElementsNearLine( const gp_Ax1&                           line,
91                                     SMDSAbs_ElementType                     type,
92                                     std::vector< const SMDS_MeshElement* >& foundElems) = 0;
93   /*!
94    * \brief Return elements whose bounding box intersects a sphere
95    */
96   virtual void GetElementsInSphere( const gp_XYZ&                           center,
97                                     const double                            radius,
98                                     SMDSAbs_ElementType                     type,
99                                     std::vector< const SMDS_MeshElement* >& foundElems) = 0;
100   /*!
101    * \brief Return elements whose bounding box intersects a given bounding box
102    */
103   virtual void GetElementsInBox( const Bnd_B3d&                          box,
104                                  SMDSAbs_ElementType                     type,
105                                  std::vector< const SMDS_MeshElement* >& foundElems) = 0;
106   /*!
107    * \brief Find out if the given point is out of closed 2D mesh.
108    */
109   virtual TopAbs_State GetPointState(const gp_Pnt& point) = 0;
110
111   /*!
112    * \brief Return a projection of a given point to a 2D mesh.
113    *        Optionally return the closest face
114    */
115   virtual gp_XYZ Project(const gp_Pnt&            point,
116                          SMDSAbs_ElementType      type,
117                          const SMDS_MeshElement** closestFace= 0) = 0;
118
119   virtual ~SMESH_ElementSearcher();
120 };
121
122 namespace SMESH_MeshAlgos
123 {
124   /*!
125    * \brief Return SMESH_NodeSearcher. The caller is responsible for deleting it
126    */
127   SMESHUtils_EXPORT
128   SMESH_NodeSearcher* GetNodeSearcher( SMDS_Mesh& mesh );
129
130   SMESHUtils_EXPORT
131   SMESH_NodeSearcher* GetNodeSearcher( SMDS_ElemIteratorPtr elemIt );
132
133   /*!
134    * \brief Return SMESH_ElementSearcher. The caller is responsible for deleting it
135    */
136   SMESHUtils_EXPORT
137   SMESH_ElementSearcher* GetElementSearcher( SMDS_Mesh& mesh,
138                                              double     tolerance=-1.);
139   SMESHUtils_EXPORT
140   SMESH_ElementSearcher* GetElementSearcher( SMDS_Mesh& mesh,
141                                              SMDS_ElemIteratorPtr elemIt,
142                                              double     tolerance=-1. );
143
144
145   /*!
146    * \brief Return true if the point is IN or ON of the element
147    */
148   SMESHUtils_EXPORT
149   bool IsOut( const SMDS_MeshElement* element, const gp_Pnt& point, double tol );
150
151   SMESHUtils_EXPORT
152   double GetDistance( const SMDS_MeshElement* elem, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
153
154   SMESHUtils_EXPORT
155   double GetDistance( const SMDS_MeshEdge* edge, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
156
157   SMESHUtils_EXPORT
158   double GetDistance( const SMDS_MeshFace* face, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
159
160   SMESHUtils_EXPORT
161   double GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt& point, gp_XYZ* closestPnt = 0 );
162
163   SMESHUtils_EXPORT
164   void GetBarycentricCoords( const gp_XY& point,
165                              const gp_XY& t0, const gp_XY& t1, const gp_XY& t2,
166                              double &    bc0, double &    bc1);
167
168   /*!
169    * Return a face having linked nodes n1 and n2 and which is
170    * - not in avoidSet,
171    * - in elemSet provided that !elemSet.empty()
172    * i1 and i2 optionally returns indices of n1 and n2
173    */
174   SMESHUtils_EXPORT
175   const SMDS_MeshElement* FindFaceInSet(const SMDS_MeshNode*    n1,
176                                         const SMDS_MeshNode*    n2,
177                                         const TIDSortedElemSet& elemSet,
178                                         const TIDSortedElemSet& avoidSet,
179                                         int*                    i1=0,
180                                         int*                    i2=0);
181   /*!
182    * \brief Calculate normal of a mesh face
183    */
184   SMESHUtils_EXPORT
185   bool FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool normalized=true);
186
187   /*!
188    * \brief Return nodes common to two elements
189    */
190   SMESHUtils_EXPORT
191   std::vector< const SMDS_MeshNode*> GetCommonNodes(const SMDS_MeshElement* e1,
192                                                     const SMDS_MeshElement* e2);
193   /*!
194    * \brief Return true if node1 encounters first in the face and node2, after.
195    *        The nodes are supposed to be neighbor nodes in the face.
196    */
197   SMESHUtils_EXPORT
198   bool IsRightOrder( const SMDS_MeshElement* face,
199                      const SMDS_MeshNode*    node0,
200                      const SMDS_MeshNode*    node1 );
201   /*!
202    * \brief Mark elements given by SMDS_Iterator
203    */
204   template< class ElemIter >
205   void MarkElems( ElemIter it, const bool isMarked )
206   {
207     while ( it->more() ) it->next()->setIsMarked( isMarked );
208   }
209   /*!
210    * \brief Mark elements given by std iterators
211    */
212   template< class ElemIter >
213   void MarkElems( ElemIter it, ElemIter end, const bool isMarked )
214   {
215     for ( ; it != end; ++it ) (*it)->setIsMarked( isMarked );
216   }
217   /*!
218    * \brief Mark nodes of elements given by SMDS_Iterator
219    */
220   template< class ElemIter >
221   void MarkElemNodes( ElemIter it, const bool isMarked, const bool markElem = false )
222   {
223     if ( markElem )
224       while ( it->more() ) {
225         const SMDS_MeshElement* e = it->next();
226         e->setIsMarked( isMarked );
227         MarkElems( e->nodesIterator(), isMarked );
228       }
229     else
230       while ( it->more() )
231         MarkElems( it->next()->nodesIterator(), isMarked );
232   }
233   /*!
234    * \brief Mark elements given by std iterators
235    */
236   template< class ElemIter >
237   void MarkElemNodes( ElemIter it, ElemIter end, const bool isMarked, const bool markElem = false )
238   {
239     if ( markElem )
240       for ( ; it != end; ++it ) {
241         (*it)->setIsMarked( isMarked );
242         MarkElems( (*it)->nodesIterator(), isMarked );
243       }
244     else
245       for ( ; it != end; ++it )
246         MarkElems( (*it)->nodesIterator(), isMarked );
247   }
248
249   // 2 nodes + optional medium node
250   struct Edge
251   {
252     const SMDS_MeshNode* _node1;
253     const SMDS_MeshNode* _node2;
254     const SMDS_MeshNode* _medium;
255   };
256
257   /*!
258    * Return sharp edges of faces and non-manifold ones.
259    * Optionally adds existing edges to the result. Angle is in degrees.
260    */
261   SMESHUtils_EXPORT
262   std::vector< Edge > FindSharpEdges( SMDS_Mesh* mesh,
263                                       double     angle,
264                                       bool       addExisting );
265
266   /*!
267    * Distribute all faces of the mesh between groups using given edges.
268    */
269   SMESHUtils_EXPORT
270   std::vector< std::vector< const SMDS_MeshElement* > >
271   SeparateFacesByEdges( SMDS_Mesh* mesh, const std::vector< Edge >& edges );
272
273
274   typedef std::vector<const SMDS_MeshNode*> TFreeBorder;
275   typedef std::vector<TFreeBorder>          TFreeBorderVec;
276   struct TFreeBorderPart
277   {
278     int _border; // border index within a TFreeBorderVec
279     int _node1;  // node index within the border-th TFreeBorder
280     int _node2;
281     int _nodeLast;
282   };
283   typedef std::vector<TFreeBorderPart>  TCoincidentGroup;
284   typedef std::vector<TCoincidentGroup> TCoincidentGroupVec;
285   struct CoincidentFreeBorders
286   {
287     TFreeBorderVec      _borders;          // nodes of all free borders
288     TCoincidentGroupVec _coincidentGroups; // groups of coincident parts of borders
289   };
290
291   /*!
292    * Returns TFreeBorder's coincident within the given tolerance.
293    * If the tolerance <= 0.0 then one tenth of an average size of elements adjacent
294    * to free borders being compared is used.
295    */
296   SMESHUtils_EXPORT
297   void FindCoincidentFreeBorders(SMDS_Mesh&              mesh,
298                                  double                  tolerance,
299                                  CoincidentFreeBorders & foundFreeBordes);
300   // Implemented in ./SMESH_FreeBorders.cxx
301
302   /*!
303    * Returns all or only closed TFreeBorder's.
304    * Optionally check if the mesh is manifold and if faces are correctly oriented.
305    */
306   SMESHUtils_EXPORT
307   void FindFreeBorders(SMDS_Mesh&       mesh,
308                        TFreeBorderVec & foundFreeBordes,
309                        const bool       closedOnly,
310                        bool*            isManifold = 0,
311                        bool*            isGoodOri = 0);
312   // Implemented in ./SMESH_FreeBorders.cxx
313
314   /*!
315    * Fill a hole defined by a TFreeBorder with 2D elements.
316    */
317   SMESHUtils_EXPORT
318   void FillHole(const TFreeBorder &                   freeBorder,
319                 SMDS_Mesh&                            mesh,
320                 std::vector<const SMDS_MeshElement*>& newFaces);
321   // Implemented in ./SMESH_FillHole.cxx
322
323
324   /*!
325    * \brief Find nodes whose merge makes the element invalid
326    */
327   SMESHUtils_EXPORT
328   void DeMerge(const SMDS_MeshElement*              elem,
329                std::vector< const SMDS_MeshNode* >& newNodes,
330                std::vector< const SMDS_MeshNode* >& noMergeNodes);
331   // Implemented in SMESH_DeMerge.cxx
332
333
334   typedef std::vector< std::pair< const SMDS_MeshElement*, const SMDS_MeshElement* > > TEPairVec;
335   typedef std::vector< std::pair< const SMDS_MeshNode*, const SMDS_MeshNode* > >       TNPairVec;
336   /*!
337    * \brief Create an offset mesh of given faces
338    *  \param [in] faceIt - the input faces
339    *  \param [in] theFixIntersections - to fix self intersections of the offset mesh or not
340    *  \param [out] new2OldFaces - history of faces
341    *  \param [out] new2OldNodes - history of nodes
342    *  \return SMDS_Mesh* - the new offset mesh, a caller should delete
343    */
344   SMESHUtils_EXPORT
345   SMDS_Mesh* MakeOffset( SMDS_ElemIteratorPtr faceIt,
346                          SMDS_Mesh&           mesh,
347                          const double         offset,
348                          const bool           theFixIntersections,
349                          TEPairVec&           new2OldFaces,
350                          TNPairVec&           new2OldNodes );
351   // Implemented in ./SMESH_Offset.cxx
352
353
354   /*!
355    * \brief Divide a mesh face into triangles
356    */
357   // Implemented in ./SMESH_Triangulate.cxx
358
359   class SMESHUtils_EXPORT Triangulate
360   {
361   public:
362
363     static int GetNbTriangles( const SMDS_MeshElement* face );
364
365     int GetTriangles( const SMDS_MeshElement*             face,
366                       std::vector< const SMDS_MeshNode*>& nodes);
367   private:
368
369     bool triangulate( std::vector< const SMDS_MeshNode*>& nodes, const size_t nbNodes );
370
371     /*!
372      * \brief Vertex of a polygon. Together with 2 neighbor Vertices represents a triangle
373      */
374     struct PolyVertex
375     {
376       SMESH_NodeXYZ _nxyz;
377       gp_XY         _xy;
378       PolyVertex*   _prev;
379       PolyVertex*   _next;
380
381       void   SetNodeAndNext( const SMDS_MeshNode* n, PolyVertex& v );
382       void   GetTriaNodes( const SMDS_MeshNode** nodes) const;
383       double TriaArea() const;
384       bool   IsInsideTria( const PolyVertex* v );
385       PolyVertex* Delete();
386     };
387     std::vector< PolyVertex > _pv;
388   };
389
390
391 } // namespace SMESH_MeshAlgos
392
393 #endif