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