Salome HOME
23491: EDF 15591 - Duplicate Elements / Nodes
[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
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 SMDS_MeshNode;
45 class SMDS_MeshElement;
46 class SMDS_Mesh;
47
48 //=======================================================================
49 /*!
50  * \brief Searcher for the node closest to a point
51  */
52 //=======================================================================
53
54 struct SMESHUtils_EXPORT SMESH_NodeSearcher
55 {
56   virtual const SMDS_MeshNode* FindClosestTo( const gp_Pnt& pnt ) = 0;
57   virtual void MoveNode( const SMDS_MeshNode* node, const gp_Pnt& toPnt ) = 0;
58   virtual int  FindNearPoint(const gp_Pnt&                        point,
59                              const double                         tolerance,
60                              std::vector< const SMDS_MeshNode* >& foundNodes) = 0;
61   virtual ~SMESH_NodeSearcher() {}
62 };
63
64 //=======================================================================
65 /*!
66  * \brief Searcher for elements
67  */
68 //=======================================================================
69
70 struct SMESHUtils_EXPORT SMESH_ElementSearcher
71 {
72   /*!
73    * \brief Find elements of given type where the given point is IN or ON.
74    *        Returns nb of found elements and elements them-selves.
75    *
76    * 'ALL' type means elements of any type excluding nodes and 0D elements
77    */
78   virtual int FindElementsByPoint(const gp_Pnt&                           point,
79                                   SMDSAbs_ElementType                     type,
80                                   std::vector< const SMDS_MeshElement* >& foundElems) = 0;
81   /*!
82    * \brief Return an element most close to the given point
83    */
84   virtual const SMDS_MeshElement* FindClosestTo( const gp_Pnt&       point,
85                                                  SMDSAbs_ElementType type) = 0;
86   /*!
87    * \brief Return elements possibly intersecting the line
88    */
89   virtual void GetElementsNearLine( const gp_Ax1&                           line,
90                                     SMDSAbs_ElementType                     type,
91                                     std::vector< const SMDS_MeshElement* >& foundElems) = 0;
92   /*!
93    * \brief Return elements whose bounding box intersects a sphere
94    */
95   virtual void GetElementsInSphere( const gp_XYZ&                           center,
96                                     const double                            radius,
97                                     SMDSAbs_ElementType                     type,
98                                     std::vector< const SMDS_MeshElement* >& foundElems) = 0;
99   /*!
100    * \brief Find out if the given point is out of closed 2D mesh.
101    */
102   virtual TopAbs_State GetPointState(const gp_Pnt& point) = 0;
103   virtual ~SMESH_ElementSearcher();
104 };
105
106 namespace SMESH_MeshAlgos
107 {
108   /*!
109    * \brief Return true if the point is IN or ON of the element
110    */
111   SMESHUtils_EXPORT
112   bool IsOut( const SMDS_MeshElement* element, const gp_Pnt& point, double tol );
113
114   SMESHUtils_EXPORT
115   double GetDistance( const SMDS_MeshElement* elem, const gp_Pnt& point );
116
117   SMESHUtils_EXPORT
118   double GetDistance( const SMDS_MeshEdge* edge, const gp_Pnt& point );
119
120   SMESHUtils_EXPORT
121   double GetDistance( const SMDS_MeshFace* face, const gp_Pnt& point );
122
123   SMESHUtils_EXPORT
124   double GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt& point );
125
126   SMESHUtils_EXPORT
127   void GetBarycentricCoords( const gp_XY& point,
128                              const gp_XY& t0, const gp_XY& t1, const gp_XY& t2,
129                              double &    bc0, double &    bc1);
130
131   /*!
132    * Return a face having linked nodes n1 and n2 and which is
133    * - not in avoidSet,
134    * - in elemSet provided that !elemSet.empty()
135    * i1 and i2 optionally returns indices of n1 and n2
136    */
137   SMESHUtils_EXPORT
138   const SMDS_MeshElement* FindFaceInSet(const SMDS_MeshNode*    n1,
139                                         const SMDS_MeshNode*    n2,
140                                         const TIDSortedElemSet& elemSet,
141                                         const TIDSortedElemSet& avoidSet,
142                                         int*                    i1=0,
143                                         int*                    i2=0);
144   /*!
145    * \brief Calculate normal of a mesh face
146    */
147   SMESHUtils_EXPORT
148   bool FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool normalized=true);
149
150   /*!
151    * \brief Return nodes common to two elements
152    */
153   SMESHUtils_EXPORT
154   std::vector< const SMDS_MeshNode*> GetCommonNodes(const SMDS_MeshElement* e1,
155                                                     const SMDS_MeshElement* e2);
156
157   /*!
158    * \brief Mark elements given by SMDS_Iterator
159    */
160   template< class ElemIter >
161   void MarkElems( ElemIter it, const bool isMarked )
162   {
163     while ( it->more() ) it->next()->setIsMarked( isMarked );
164   }
165   /*!
166    * \brief Mark elements given by std iterators
167    */
168   template< class ElemIter >
169   void MarkElems( ElemIter it, ElemIter end, const bool isMarked )
170   {
171     for ( ; it != end; ++it ) (*it)->setIsMarked( isMarked );
172   }
173   /*!
174    * \brief Mark nodes of elements given by SMDS_Iterator
175    */
176   template< class ElemIter >
177   void MarkElemNodes( ElemIter it, const bool isMarked, const bool markElem = false )
178   {
179     if ( markElem )
180       while ( it->more() ) {
181         const SMDS_MeshElement* e = it->next();
182         e->setIsMarked( isMarked );
183         MarkElems( e->nodesIterator(), isMarked );
184       }
185     else
186       while ( it->more() )
187         MarkElems( it->next()->nodesIterator(), isMarked );
188   }
189   /*!
190    * \brief Mark elements given by std iterators
191    */
192   template< class ElemIter >
193   void MarkElemNodes( ElemIter it, ElemIter end, const bool isMarked, const bool markElem = false )
194   {
195     if ( markElem )
196       for ( ; it != end; ++it ) {
197         (*it)->setIsMarked( isMarked );
198         MarkElems( (*it)->nodesIterator(), isMarked );
199       }
200     else
201       for ( ; it != end; ++it )
202         MarkElems( (*it)->nodesIterator(), isMarked );
203   }
204
205   /*!
206    * \brief Return SMESH_NodeSearcher. The caller is responsible for deleteing it
207    */
208   SMESHUtils_EXPORT
209   SMESH_NodeSearcher* GetNodeSearcher( SMDS_Mesh& mesh );
210
211   SMESHUtils_EXPORT
212   SMESH_NodeSearcher* GetNodeSearcher( SMDS_ElemIteratorPtr elemIt );
213
214   /*!
215    * \brief Return SMESH_ElementSearcher. The caller is responsible for deleting it
216    */
217   SMESHUtils_EXPORT
218   SMESH_ElementSearcher* GetElementSearcher( SMDS_Mesh& mesh,
219                                              double     tolerance=-1.);
220   SMESHUtils_EXPORT
221   SMESH_ElementSearcher* GetElementSearcher( SMDS_Mesh& mesh,
222                                              SMDS_ElemIteratorPtr elemIt,
223                                              double     tolerance=-1. );
224
225
226
227   typedef std::vector<const SMDS_MeshNode*> TFreeBorder;
228   typedef std::vector<TFreeBorder>          TFreeBorderVec;
229   struct TFreeBorderPart
230   {
231     int _border; // border index within a TFreeBorderVec
232     int _node1;  // node index within the border-th TFreeBorder
233     int _node2;
234     int _nodeLast;
235   };
236   typedef std::vector<TFreeBorderPart>  TCoincidentGroup;
237   typedef std::vector<TCoincidentGroup> TCoincidentGroupVec;
238   struct CoincidentFreeBorders
239   {
240     TFreeBorderVec      _borders;          // nodes of all free borders
241     TCoincidentGroupVec _coincidentGroups; // groups of coincident parts of borders
242   };
243
244   /*!
245    * Returns TFreeBorder's coincident within the given tolerance.
246    * If the tolerance <= 0.0 then one tenth of an average size of elements adjacent
247    * to free borders being compared is used.
248    *
249    * (Implemented in ./SMESH_FreeBorders.cxx)
250    */
251   SMESHUtils_EXPORT
252   void FindCoincidentFreeBorders(SMDS_Mesh&              mesh,
253                                  double                  tolerance,
254                                  CoincidentFreeBorders & foundFreeBordes);
255   
256
257   /*!
258    * \brief Find nodes whose merge makes the element invalid
259    *
260    * (Implemented in SMESH_DeMerge.cxx)
261    */
262   SMESHUtils_EXPORT
263   void DeMerge(const SMDS_MeshElement*              elem,
264                std::vector< const SMDS_MeshNode* >& newNodes,
265                std::vector< const SMDS_MeshNode* >& noMergeNodes);
266
267 } // namespace SMESH_MeshAlgos
268
269 #endif