Salome HOME
aa38284431b6bb048aa943f5aeeef368792cd499
[modules/smesh.git] / src / StdMeshers / StdMeshers_ProjectionUtils.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
23 //  SMESH SMESH : idl implementation based on 'SMESH' unit's calsses
24 // File      : StdMeshers_ProjectionUtils.hxx
25 // Created   : Thu Oct 26 15:37:24 2006
26 // Author    : Edward AGAPOV (eap)
27 //
28 #ifndef StdMeshers_ProjectionUtils_HeaderFile
29 #define StdMeshers_ProjectionUtils_HeaderFile
30
31 #include "SMESH_StdMeshers.hxx"
32
33 #include "SMDS_MeshElement.hxx"
34 #include "SMESH_Delaunay.hxx"
35 #include "StdMeshers_FaceSide.hxx"
36
37 #include <ShapeAnalysis_Surface.hxx>
38 #include <TopTools_DataMapOfShapeShape.hxx>
39 #include <TopTools_IndexedDataMapOfShapeListOfShape.hxx>
40 #include <TopTools_IndexedMapOfShape.hxx>
41 #include <TopoDS_Edge.hxx>
42 #include <TopoDS_Face.hxx>
43 #include <TopoDS_Vertex.hxx>
44 #include <gp_GTrsf.hxx>
45 #include <gp_GTrsf2d.hxx>
46
47 #include <list>
48 #include <map>
49
50 class SMDS_MeshNode;
51 class SMESH_Algo;
52 class SMESH_Hypothesis;
53 class SMESH_Mesh;
54 class SMESH_subMesh;
55 class TopoDS_Shape;
56
57 //-----------------------------------------------------------------------------------------
58 /*!
59  * \brief Struct used instead of a sole TopTools_DataMapOfShapeShape to avoid
60  *        problems with bidirectional bindings
61  */
62 struct StdMeshers_ShapeShapeBiDirectionMap
63 {
64   TopTools_DataMapOfShapeShape _map1to2, _map2to1;
65
66   enum EAssocType {
67     UNDEF, INIT_VERTEX, PROPAGATION, PARTNER, CLOSE_VERTEX, COMMON_VERTEX, FEW_EF };
68   EAssocType _assocType;
69
70   // convention: s1 - target, s2 - source
71   bool Bind( const TopoDS_Shape& s1, const TopoDS_Shape& s2 )
72   { _map1to2.Bind( s1, s2 ); return _map2to1.Bind( s2, s1 ); }
73   bool IsBound( const TopoDS_Shape& s, const bool isShape2=false ) const 
74   { return (isShape2 ? _map2to1 : _map1to2).IsBound( s ); }
75   bool IsEmpty() const { return _map1to2.IsEmpty(); }
76   int  Extent()  const { return _map1to2.Extent(); }
77   void Clear() { _map1to2.Clear(); _map2to1.Clear(); }
78   const TopoDS_Shape& operator()( const TopoDS_Shape& s, const bool isShape2=false ) const
79   { // if we get a Standard_NoSuchObject here, it means that the calling code
80     // passes incorrect isShape2
81     return (isShape2 ? _map2to1 : _map1to2)( s );
82   }
83   StdMeshers_ShapeShapeBiDirectionMap() : _assocType( UNDEF ) {}
84   void SetAssocType( EAssocType type ) { if ( _assocType == UNDEF ) _assocType = type; }
85 };
86
87 /*!
88  * \brief Methods common to Projection algorithms
89  */
90 namespace StdMeshers_ProjectionUtils
91 {
92   typedef StdMeshers_ShapeShapeBiDirectionMap                  TShapeShapeMap;
93   typedef TopTools_IndexedDataMapOfShapeListOfShape            TAncestorMap;
94   typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*,
95                    TIDCompare>                                 TNodeNodeMap;
96
97
98   //-----------------------------------------------------------------------------------------
99   /*!
100    * \brief Finds transformation between two sets of 2D points using
101    *        a least square approximation
102    */
103   class TrsfFinder2D
104   {
105     gp_GTrsf2d _trsf;
106     gp_XY      _srcOrig;
107   public:
108     TrsfFinder2D(): _srcOrig(0,0) {}
109
110     void Set( const gp_GTrsf2d& t ) { _trsf = t; } // it's an alternative to Solve()
111
112     bool Solve( const std::vector< gp_XY >& srcPnts,
113                 const std::vector< gp_XY >& tgtPnts );
114
115     gp_XY Transform( const gp_Pnt2d& srcUV ) const;
116
117     bool IsIdentity() const { return ( _trsf.Form() == gp_Identity ); }
118   };
119   //-----------------------------------------------------------------------------------------
120   /*!
121    * \brief Finds transformation between two sets of 3D points using
122    *        a least square approximation
123    */
124   class TrsfFinder3D
125   {
126     gp_GTrsf _trsf;
127     gp_XYZ   _srcOrig;
128   public:
129     TrsfFinder3D(): _srcOrig(0,0,0) {}
130
131     void Set( const gp_GTrsf& t ) { _trsf = t; } // it's an alternative to Solve()
132
133     bool Solve( const std::vector< gp_XYZ > & srcPnts,
134                 const std::vector< gp_XYZ > & tgtPnts );
135
136     gp_XYZ Transform( const gp_Pnt& srcP ) const;
137
138     gp_XYZ TransformVec( const gp_Vec& v ) const;
139
140     bool IsIdentity() const { return ( _trsf.Form() == gp_Identity ); }
141
142     bool Invert();
143   };
144
145   //-----------------------------------------------------------------------------------------
146   /*!
147    * \brief Create a Delaunay triangulation of nodes on a face boundary
148    *        and provide exploration of nodes shared by elements lying on
149    *        the face. For a returned node, also return a Delaunay triangle
150    *        the node lies in and its Barycentric Coordinates within the triangle.
151    *        Only non-marked nodes are visited.
152    *
153    * The main methods are defined in ../SMESHUtils/SMESH_Delaunay.hxx
154    */
155   class Delaunay : public SMESH_Delaunay
156   {
157   public:
158
159     Delaunay( const TSideVector& wires, bool checkUV = false );
160
161     Delaunay( const std::vector< const UVPtStructVec* > & boundaryNodes,
162               SMESH_MesherHelper&                         faceHelper,
163               bool                                        checkUV = false);
164
165   protected:
166     virtual gp_XY getNodeUV( const TopoDS_Face& face, const SMDS_MeshNode* node ) const;
167
168   private:
169     SMESH_MesherHelper*    _helper;
170     StdMeshers_FaceSidePtr _wire;
171     bool                  *_checkUVPtr, _checkUV;
172   };
173   typedef boost::shared_ptr< Delaunay > DelaunayPtr;
174
175   //-----------------------------------------------------------------------------------------
176   /*!
177    * \brief Morph mesh on the target FACE to lie within FACE boundary w/o distortion
178    */
179   class Morph
180   {
181     Delaunay       _delaunay;
182     SMESH_subMesh* _srcSubMesh;
183   public:
184
185     Morph(const TSideVector& srcWires);
186
187     bool Perform(SMESH_MesherHelper&           tgtHelper,
188                  const TSideVector&            tgtWires,
189                  Handle(ShapeAnalysis_Surface) tgtSurface,
190                  const TNodeNodeMap&           src2tgtNodes,
191                  const bool                    moveAll);
192   };
193
194   //-----------------------------------------------------------------------------------------
195   /*!
196    * \brief Looks for association of all sub-shapes of two shapes
197    * \param theShape1 - shape 1
198    * \param theMesh1 - mesh built on shape 1
199    * \param theShape2 - shape 2
200    * \param theMesh2 - mesh built on shape 2
201    * \param theAssociation - association map to be filled that may
202    *                         contain association of one or two pairs of vertices
203    * \retval bool - true if association found
204    */
205   bool FindSubShapeAssociation(const TopoDS_Shape& theShape1,
206                                SMESH_Mesh*         theMesh1,
207                                const TopoDS_Shape& theShape2,
208                                SMESH_Mesh*         theMesh2,
209                                TShapeShapeMap &    theAssociationMap);
210
211   /*!
212    * \brief Find association of edges of faces
213    *  \param face1 - face 1
214    *  \param VV1 - vertices of face 1
215    *  \param face2 - face 2
216    *  \param VV2 - vertices of face 2 associated with oned of face 1
217    *  \param edges1 - out list of edges of face 1
218    *  \param edges2 - out list of edges of face 2
219    *  \param isClosenessAssoc - is association starting by VERTEX closeness
220    *  \retval int - nb of edges in an outer wire in a success case, else zero
221    */
222   int FindFaceAssociation(const TopoDS_Face&         face1,
223                           TopoDS_Vertex              VV1[2],
224                           const TopoDS_Face&         face2,
225                           TopoDS_Vertex              VV2[2],
226                           std::list< TopoDS_Edge > & edges1,
227                           std::list< TopoDS_Edge > & edges2,
228                           const bool                 isClosenessAssoc=false);
229
230   /*!
231    * \brief Insert vertex association defined by a hypothesis into a map
232    * \param theHyp - hypothesis
233    * \param theAssociationMap - association map
234    * \param theTargetShape - the shape theHyp assigned to
235    */
236   void InitVertexAssociation( const SMESH_Hypothesis* theHyp,
237                               TShapeShapeMap &        theAssociationMap);
238
239   /*!
240    * \brief Inserts association theShape1 <-> theShape2 to TShapeShapeMap
241    * \param theShape1 - target shape
242    * \param theShape2 - source shape
243    * \param theAssociationMap - association map 
244    * \param theBidirectional - if false, inserts theShape1 -> theShape2 association
245    * \retval bool - true if there was no association for these shapes before
246    */
247   bool InsertAssociation( const TopoDS_Shape& theShape1, // target
248                           const TopoDS_Shape& theShape2, // source
249                           TShapeShapeMap &    theAssociationMap);
250
251   /*!
252    * \brief Finds an edge by its vertices in a main shape of the mesh
253    */
254   TopoDS_Edge GetEdgeByVertices( SMESH_Mesh*          aMesh,
255                                  const TopoDS_Vertex& V1,
256                                  const TopoDS_Vertex& V2);
257
258   /*!
259    * \brief Return another face sharing an edge
260    * \param edgeToFaces - data map of descendants to ancestors
261    */
262   TopoDS_Face GetNextFace( const TAncestorMap& edgeToFaces,
263                            const TopoDS_Edge&  edge,
264                            const TopoDS_Face&  face);
265   /*!
266    * \brief Return other vertex of an edge
267    */
268   TopoDS_Vertex GetNextVertex(const TopoDS_Edge&   edge,
269                               const TopoDS_Vertex& vertex);
270
271   /*!
272    * \brief Return an oriented propagation edge
273    * \param aMesh - mesh
274    * \param fromEdge - start edge for propagation
275    * \param chain - return, if provided, a propagation chain passed till
276    *        anEdge; if anEdge.IsNull() then a full propagation chain is returned
277    * \retval pair<int,TopoDS_Edge> - propagation step and found edge
278    */
279   std::pair<int,TopoDS_Edge> GetPropagationEdge( SMESH_Mesh*                 aMesh,
280                                                  const TopoDS_Edge&          anEdge,
281                                                  const TopoDS_Edge&          fromEdge,
282                                                  TopTools_IndexedMapOfShape* chain=0);
283
284   /*!
285    * \brief Find corresponding nodes on two faces
286    * \param face1 - the first face
287    * \param mesh1 - mesh containing elements on the first face
288    * \param face2 - the second face
289    * \param mesh2 - mesh containing elements on the second face
290    * \param assocMap - map associating sub-shapes of the faces
291    * \param nodeIn2OutMap - map containing found matching nodes
292    * \retval bool - is a success
293    */
294   bool FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
295                                  SMESH_Mesh*            mesh1,
296                                  const TopoDS_Face&     face2,
297                                  SMESH_Mesh*            mesh2,
298                                  const TShapeShapeMap & assocMap,
299                                  TNodeNodeMap &         nodeIn2OutMap);
300   /*!
301    * \brief Return any sub-shape of a face belonging to the outer wire
302    * \param face - the face
303    * \param type - type of sub-shape to return
304    * \retval TopoDS_Shape - the found sub-shape
305    */
306   TopoDS_Shape OuterShape( const TopoDS_Face& face,
307                            TopAbs_ShapeEnum   type);
308
309   /*!
310    * \brief Check that submeshis is computed and try to compute it if is not
311    * \param sm - submesh to compute
312    * \param iterationNb - int used to stop infinite recursive call
313    * \retval bool - true if computed
314    */
315   bool MakeComputed(SMESH_subMesh * sm, const int iterationNb = 0);
316
317   /*!
318    * \brief Returns an error message to show in case if MakeComputed( sm ) fails.
319    */
320   std::string SourceNotComputedError( SMESH_subMesh * sm = 0,
321                                       SMESH_Algo*     projAlgo=0);
322
323   /*!
324    * \brief Set event listeners to submesh with projection algo
325    * \param subMesh - submesh with projection algo
326    * \param srcShape - source shape
327    * \param srcMesh - source mesh
328    */
329   void SetEventListener(SMESH_subMesh* subMesh,
330                         TopoDS_Shape   srcShape,
331                         SMESH_Mesh*    srcMesh);
332
333   /*!
334    * \brief Return a boundary EDGE (or all boundary EDGEs) of edgeContainer
335    */
336   TopoDS_Edge GetBoundaryEdge(const TopoDS_Shape&       edgeContainer,
337                               const SMESH_Mesh&         mesh,
338                               std::list< TopoDS_Edge >* allBndEdges = 0 );
339 };
340
341 #endif