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