Salome HOME
Merge from V6_3_BR 06/06/2011
[modules/smesh.git] / src / StdMeshers / StdMeshers_Prism_3D.hxx
1 // Copyright (C) 2007-2011  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.
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 : implementaion of SMESH idl descriptions
24 //  File   : StdMeshers_Prism_3D.hxx
25 //  Module : SMESH
26 //
27 #ifndef _SMESH_Prism_3D_HXX_
28 #define _SMESH_Prism_3D_HXX_
29
30 #include "SMESH_StdMeshers.hxx"
31
32 #include "SMESH_3D_Algo.hxx"
33 #include "SMDS_TypeOfPosition.hxx"
34 #include "SMDS_MeshNode.hxx"
35 #include "SMESH_Block.hxx"
36 #include "SMESH_Mesh.hxx"
37 #include "SMESHDS_Mesh.hxx"
38 #include "SMESH_subMesh.hxx"
39 #include "SMESH_MesherHelper.hxx"
40 #include "SMESH_Comment.hxx"
41
42 #include <vector>
43
44 #include <Adaptor3d_Curve.hxx>
45 #include <Adaptor3d_Surface.hxx>
46 #include <Adaptor2d_Curve2d.hxx>
47 #include <BRepAdaptor_Surface.hxx>
48 #include <TopTools_IndexedMapOfOrientedShape.hxx>
49 #include <gp_XYZ.hxx>
50 #include <gp_Trsf.hxx>
51
52
53 class SMESHDS_SubMesh;
54 class TopoDS_Edge;
55 class TopoDS_Faces;
56 struct TNode;
57
58 //typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*> TNodeNodeMap;
59 typedef std::vector<const SMDS_MeshNode* > TNodeColumn;
60
61 // map of bottom nodes to the column of nodes above them
62 // (the column includes the bottom nodes)
63 typedef std::map< TNode, TNodeColumn >  TNode2ColumnMap;
64 typedef std::map< double, TNodeColumn > TParam2ColumnMap;
65 typedef std::map< double, TNodeColumn >::const_iterator TParam2ColumnIt;
66
67 typedef TopTools_IndexedMapOfOrientedShape TBlockShapes;
68
69 // ===============================================
70 /*!
71  * \brief Structure containing node relative data
72  */
73 // ===============================================
74
75 struct TNode
76 {
77   const SMDS_MeshNode* myNode;
78   mutable gp_XYZ       myParams;
79
80   gp_XYZ GetCoords() const { return gp_XYZ( myNode->X(), myNode->Y(), myNode->Z() ); }
81   gp_XYZ GetParams() const { return myParams; }
82   gp_XYZ& ChangeParams() const { return myParams; }
83   bool HasParams() const { return myParams.X() >= 0.0; }
84   SMDS_TypeOfPosition GetPositionType() const
85   { return myNode ? myNode->GetPosition()->GetTypeOfPosition() : SMDS_TOP_UNSPEC; }
86   bool IsNeighbor( const TNode& other ) const;
87
88   TNode(const SMDS_MeshNode* node = 0): myNode(node), myParams(-1,-1,-1) {}
89   bool operator < (const TNode& other) const { return myNode->GetID() < other.myNode->GetID(); }
90 };
91
92 // ===============================================================
93 /*!
94  * \brief Tool analyzing and giving access to a prism geometry 
95  *  treating it like a block, i.e. the four side faces are
96  *  emulated by division/uniting of missing/excess faces.
97  *  It also manage associations between block subshapes and a mesh.
98  */
99 // ===============================================================
100
101 class STDMESHERS_EXPORT StdMeshers_PrismAsBlock: public SMESH_Block
102 {
103 public:
104   /*!
105    * \brief Constructor. Initialization is needed
106    */
107   StdMeshers_PrismAsBlock();
108
109   ~StdMeshers_PrismAsBlock();
110
111   /*!
112    * \brief Initialization.
113     * \param helper - helper loaded with mesh and 3D shape
114     * \param shape3D - a closed shell or solid
115     * \retval bool - false if a mesh or a shape are KO
116     *
117     * Analyse shape geometry and mesh.
118     * If there are triangles on one of faces, it becomes 'bottom'
119    */
120   bool Init(SMESH_MesherHelper* helper, const TopoDS_Shape& shape3D);
121
122   /*!
123    * \brief Return problem description
124    */
125   SMESH_ComputeErrorPtr GetError() const { return myError; }
126
127   /*!
128    * \brief Free allocated memory
129    */
130   void Clear();
131
132   /*!
133    * \brief Return number of nodes on every vertical edge
134     * \retval int - number of nodes including end nodes
135    */
136   int VerticalSize() const { return myParam2ColumnMaps[0].begin()->second.size(); }
137
138   bool HasNotQuadElemOnTop() const { return myNotQuadOnTop; }
139
140   /*!
141    * \brief Return pointer to column of nodes
142     * \param node - bottom node from which the returned column goes up
143     * \retval const TNodeColumn* - the found column
144    */
145   const TNodeColumn* GetNodeColumn(const SMDS_MeshNode* node) const;
146
147   /*!
148    * \brief Return TParam2ColumnMap for a base edge
149     * \param baseEdgeID - base edge SMESHDS Index
150     * \param isReverse - columns in-block orientation
151     * \retval const TParam2ColumnMap& - map
152    */
153   const TParam2ColumnMap& GetParam2ColumnMap(const int baseEdgeID,
154                                              bool &    isReverse) const
155   {
156     std::pair< TParam2ColumnMap*, bool > col_frw =
157       myShapeIndex2ColumnMap.find( baseEdgeID )->second;
158     isReverse = !col_frw.second;
159     return * col_frw.first;
160   }
161
162   /*!
163    * \brief Return transformations to get coordinates of nodes of each internal layer
164    *        by nodes of the bottom. Layer is a set of nodes at a certain step
165    *        from bottom to top.
166    */
167   bool GetLayersTransformation(std::vector<gp_Trsf> & trsf) const;
168   
169   /*!
170    * \brief Return pointer to mesh
171     * \retval SMESH_Mesh - mesh
172    */
173   SMESH_Mesh* Mesh() const { return myHelper->GetMesh(); }
174
175   /*!
176    * \brief Return pointer to mesh DS
177     * \retval SMESHDS_Mesh - mesh DS
178    */
179   SMESHDS_Mesh* MeshDS() const { return Mesh()->GetMeshDS(); }
180
181   /*!
182    * \brief Return submesh of a shape
183     * \param shapeID - shape given by in-block index
184     * \retval SMESH_subMesh* - found submesh
185    */
186   SMESH_subMesh* SubMesh(const int shapeID) const
187   { return Mesh()->GetSubMesh( Shape( shapeID )); }
188
189   /*!
190    * \brief Return submesh DS of a shape
191     * \param shapeID - shape given by in-block index
192     * \retval SMESHDS_SubMesh* - found submesh DS
193    */
194   SMESHDS_SubMesh* SubMeshDS(const int shapeID) const
195   { return SubMesh(shapeID)->GetSubMeshDS(); }
196
197   /*!
198    * \brief Return a in-block shape
199     * \param shapeID - shape given by in-block index
200     * \retval SMESHDS_SubMesh* - found submesh
201    */
202   const TopoDS_Shape& Shape(const int shapeID) const
203   { return myShapeIDMap( shapeID ); }
204
205   /*!
206    * \brief Return in-block ID of a shape
207     * \param shape - block subshape
208     * \retval int - ID or zero if the shape has no ID
209    */
210   int ShapeID(const TopoDS_Shape& shape) const
211   { return myShapeIDMap.FindIndex( shape ); }
212
213   /*!
214    * \brief Check curve orientation of a bootom edge
215    *  \param meshDS - mesh DS
216    *  \param columnsMap - node columns map of side face
217    *  \param bottomEdge - the bootom edge
218    *  \param sideFaceID - side face in-block ID
219    *  \retval bool - true if orienation coinside with in-block froward orienation
220    */
221   static bool IsForwardEdge(SMESHDS_Mesh*           meshDS,
222                             const TParam2ColumnMap& columnsMap,
223                             const TopoDS_Edge &     bottomEdge,
224                             const int               sideFaceID);
225   /*!
226    * \brief Find wall faces by bottom edges
227     * \param mesh - the mesh
228     * \param mainShape - the prism
229     * \param bottomFace - the bottom face
230     * \param bottomEdges - edges bounding the bottom face
231     * \param wallFaces - faces list to fill in
232    */
233   bool GetWallFaces( SMESH_Mesh*               mesh,
234                      const TopoDS_Shape &      mainShape,
235                      const TopoDS_Shape &      bottomFace,
236                      std::list< TopoDS_Edge >& bottomEdges,
237                      std::list< int > &        nbEInW,
238                      std::list< TopoDS_Face >& wallFaces);
239
240 private:
241
242   // --------------------------------------------------------------------
243   /*!
244    * \brief Class representing a part of a geom face or
245    * a union of seleral faces. Or just an ordinary geom face
246    *
247    * It's parametrization is within [0,1] range.
248    * It redefines Adaptor3d_Surface::Value(U,V) where U and V are within [0,1]
249    */
250   // --------------------------------------------------------------------
251   class TSideFace: public Adaptor3d_Surface
252   {
253     int                             myID; //!< in-block ID
254     // map used to find out real UV by it's normalized UV
255     TParam2ColumnMap*               myParamToColumnMap;
256     BRepAdaptor_Surface             mySurface;
257     TopoDS_Edge                     myBaseEdge;
258     // first and last normalized params and orientaion for each component or it-self
259     std::vector< std::pair< double, double> > myParams;
260     bool                            myIsForward;
261     std::vector< TSideFace* >       myComponents;
262     SMESH_MesherHelper *            myHelper;
263   public:
264     TSideFace( SMESH_MesherHelper* helper,
265                const int           faceID,
266                const TopoDS_Face&  face,
267                const TopoDS_Edge&  baseEdge,
268                TParam2ColumnMap*   columnsMap,
269                const double        first = 0.0,
270                const double        last = 1.0);
271     TSideFace( const std::vector< TSideFace* >&             components,
272                const std::vector< std::pair< double, double> > & params);
273     TSideFace( const TSideFace& other );
274     ~TSideFace();
275     bool IsComplex() const
276     { return ( NbComponents() > 0 || myParams[0].first != 0. || myParams[0].second != 1. ); }
277     int FaceID() const { return myID; }
278     TParam2ColumnMap* GetColumns() const { return myParamToColumnMap; }
279     gp_XY GetNodeUV(const TopoDS_Face& F, const SMDS_MeshNode* n) const
280     { return myHelper->GetNodeUV( F, n ); }
281     const TopoDS_Edge & BaseEdge() const { return myBaseEdge; }
282     int ColumnHeight() const {
283       if ( NbComponents() ) return GetComponent(0)->GetColumns()->begin()->second.size();
284       else                  return GetColumns()->begin()->second.size(); }
285     double GetColumns(const double U, TParam2ColumnIt & col1, TParam2ColumnIt& col2 ) const;
286     int NbComponents() const { return myComponents.size(); }
287     TSideFace* GetComponent(const int i) const { return myComponents.at( i ); }
288     void SetComponent(const int i, TSideFace* c)
289     { if ( myComponents[i] ) delete myComponents[i]; myComponents[i]=c; }
290     TSideFace* GetComponent(const double U, double& localU) const;
291     bool IsForward() const { return myIsForward; }
292     // boundary geometry for a face
293     Adaptor3d_Surface* Surface() const { return new TSideFace( *this ); }
294     bool GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const;
295     Adaptor2d_Curve2d* HorizPCurve(const bool isTop, const TopoDS_Face& horFace) const;
296     Adaptor3d_Curve* HorizCurve(const bool isTop) const;
297     Adaptor3d_Curve* VertiCurve(const bool isMax) const;
298     TopoDS_Edge GetEdge( const int edge ) const;
299     int InsertSubShapes( TBlockShapes& shapeMap ) const;
300     // redefine Adaptor methods
301     gp_Pnt Value(const Standard_Real U,const Standard_Real V) const;
302     // debug
303     void dumpNodes(int nbNodes) const;
304   };
305
306   // --------------------------------------------------------------------
307   /*!
308    * \brief Class emulating geometry of a vertical edge
309    */
310   // --------------------------------------------------------------------
311   class STDMESHERS_EXPORT TVerticalEdgeAdaptor: public Adaptor3d_Curve
312   {
313     const TNodeColumn* myNodeColumn;
314   public:
315     TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter );
316     gp_Pnt Value(const Standard_Real U) const;
317     Standard_Real FirstParameter() const { return 0; }
318     Standard_Real LastParameter() const { return 1; }
319     // debug
320     void dumpNodes(int nbNodes) const;
321   };
322
323   // --------------------------------------------------------------------
324   /*!
325    * \brief Class emulating geometry of a hirizontal edge
326    */
327   // --------------------------------------------------------------------
328   class STDMESHERS_EXPORT THorizontalEdgeAdaptor: public Adaptor3d_Curve
329   {
330     const TSideFace* mySide;
331     double           myV;
332   public:
333     THorizontalEdgeAdaptor( const TSideFace* sideFace, const bool isTop)
334       :mySide(sideFace), myV( isTop ? 1.0 : 0.0 ) {}
335     gp_Pnt Value(const Standard_Real U) const;
336     Standard_Real FirstParameter() const { return 0; }
337     Standard_Real LastParameter() const { return 1; }
338     // debug
339     void dumpNodes(int nbNodes) const;
340   };
341
342   // --------------------------------------------------------------------
343   /*!
344    * \brief Class emulating pcurve on a hirizontal face
345    */
346   // --------------------------------------------------------------------
347   class STDMESHERS_EXPORT TPCurveOnHorFaceAdaptor: public Adaptor2d_Curve2d
348   {
349     const TSideFace*  mySide;
350     int               myZ;
351     TopoDS_Face       myFace;
352   public:
353     TPCurveOnHorFaceAdaptor( const TSideFace*   sideFace,
354                              const bool         isTop,
355                              const TopoDS_Face& horFace)
356       : mySide(sideFace), myFace(horFace), myZ(isTop ? mySide->ColumnHeight() - 1 : 0 ) {}
357     gp_Pnt2d Value(const Standard_Real U) const;
358     Standard_Real FirstParameter() const { return 0; }
359     Standard_Real LastParameter() const { return 1; }
360   };
361
362   bool                  myNotQuadOnTop;
363   SMESH_MesherHelper*   myHelper;
364   TBlockShapes          myShapeIDMap;
365   SMESH_ComputeErrorPtr myError;
366
367   // container of 4 side faces
368   TSideFace*            mySide; 
369   // node columns for each base edge
370   std::vector< TParam2ColumnMap >                       myParam2ColumnMaps;
371   // to find a column for a node by edge SMESHDS Index
372   std::map< int, std::pair< TParam2ColumnMap*, bool > > myShapeIndex2ColumnMap;
373
374   /*!
375    * \brief store error and comment and then return ( error == COMPERR_OK )
376    */
377   bool error(int error, const SMESH_Comment& comment = "") {
378     myError = SMESH_ComputeError::New(error,comment);
379     return myError->IsOK();
380   }
381 };
382
383 // =============================================
384 /*!
385  * \brief Algo building prisms on a prism shape
386  */
387 // =============================================
388
389 class STDMESHERS_EXPORT StdMeshers_Prism_3D: public SMESH_3D_Algo
390 {
391 public:
392   StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen);
393   virtual ~StdMeshers_Prism_3D();
394
395   virtual bool CheckHypothesis(SMESH_Mesh&                          aMesh,
396                                const TopoDS_Shape&                  aShape,
397                                SMESH_Hypothesis::Hypothesis_Status& aStatus);
398
399   virtual bool Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape);
400
401   virtual bool Evaluate(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape,
402                         MapShapeNbElems& aResMap);
403
404   /*!
405    * \brief Enable removal of quadrangles from the bottom face and
406    * triangles creation there by projection from the top
407    * (sole face meshed with triangles is considered to be a bottom one).
408    * If there are two faces with triangles, triangles must
409    * be of the same topology, else the algo fails.
410    * The method must be called before Compute()
411    */
412   void ProjectTriangles() { myProjectTriangles = true; }
413
414   /*!
415    * \brief Create prisms
416     * \param nodeColumns - columns of nodes generated from nodes of a mesh face
417     * \param helper - helper initialized by mesh and shape to add prisms to
418    */
419   static void AddPrisms( std::vector<const TNodeColumn*> & nodeColumns,
420                          SMESH_MesherHelper*          helper);
421
422 private:
423
424   /*!
425    * \brief Find correspondence between bottom and top nodes.
426    *  If elements on the bottom and top faces are topologically different,
427    *  and projection is possible and allowed, perform the projection
428     * \retval bool - is a success or not
429    */
430   bool assocOrProjBottom2Top();
431
432   /*!
433    * \brief Remove quadrangles from the top face and
434    * create triangles there by projection from the bottom
435     * \retval bool - a success or not
436    */
437   bool projectBottomToTop();
438
439   /*!
440    * \brief Set projection coordinates of a node to a face and it's subshapes
441     * \param faceID - the face given by in-block ID
442     * \param params - node normalized parameters
443     * \retval bool - is a success
444    */
445   bool setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z );
446
447 private:
448
449   bool myProjectTriangles;
450
451   StdMeshers_PrismAsBlock myBlock;
452   SMESH_MesherHelper*     myHelper;
453
454   std::vector<gp_XYZ>     myShapeXYZ; // point on each sub-shape of the block
455
456   // map of bottom nodes to the column of nodes above them
457   // (the column includes the bottom node)
458   TNode2ColumnMap         myBotToColumnMap;
459 };
460
461 #endif