]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_Prism_3D.hxx
Salome HOME
Minor DOC changes
[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 transformations to get coordinates of nodes of each internal layer
190    *        by nodes of the bottom. Layer is a set of nodes at a certain step
191    *        from bottom to top.
192    */
193   bool GetLayersTransformation(std::vector<gp_Trsf> &      trsf,
194                                const Prism_3D::TPrismTopo& prism) const;
195   
196   /*!
197    * \brief Return pointer to mesh
198     * \retval SMESH_Mesh - mesh
199    */
200   SMESH_Mesh* Mesh() const { return myHelper->GetMesh(); }
201
202   /*!
203    * \brief Return pointer to mesh DS
204     * \retval SMESHDS_Mesh - mesh DS
205    */
206   SMESHDS_Mesh* MeshDS() const { return Mesh()->GetMeshDS(); }
207
208   /*!
209    * \brief Return submesh of a shape
210     * \param shapeID - shape given by in-block index
211     * \retval SMESH_subMesh* - found submesh
212    */
213   SMESH_subMesh* SubMesh(const int shapeID) const
214   { return Mesh()->GetSubMesh( Shape( shapeID )); }
215
216   /*!
217    * \brief Return submesh DS of a shape
218     * \param shapeID - shape given by in-block index
219     * \retval SMESHDS_SubMesh* - found submesh DS
220    */
221   SMESHDS_SubMesh* SubMeshDS(const int shapeID) const
222   { return SubMesh(shapeID)->GetSubMeshDS(); }
223
224   /*!
225    * \brief Return a in-block shape
226     * \param shapeID - shape given by in-block index
227     * \retval SMESHDS_SubMesh* - found submesh
228    */
229   const TopoDS_Shape& Shape(const int shapeID) const
230   { return myShapeIDMap( shapeID ); }
231
232   /*!
233    * \brief Return in-block ID of a shape
234     * \param shape - block sub-shape
235     * \retval int - ID or zero if the shape has no ID
236    */
237   int ShapeID(const TopoDS_Shape& shape) const
238   { return myShapeIDMap.FindIndex( shape ); }
239
240   /*!
241    * \brief Check curve orientation of a bootom edge
242    *  \param meshDS - mesh DS
243    *  \param columnsMap - node columns map of side face
244    *  \param bottomEdge - the bootom edge
245    *  \param sideFaceID - side face in-block ID
246    *  \retval bool - true if orienation coinside with in-block froward orienation
247    */
248   static bool IsForwardEdge(SMESHDS_Mesh*           meshDS,
249                             const TParam2ColumnMap& columnsMap,
250                             const TopoDS_Edge &     bottomEdge,
251                             const int               sideFaceID);
252
253 private:
254
255   // --------------------------------------------------------------------
256   /*!
257    * \brief Class representing a part of a geom face or
258    * a union of seleral faces. Or just an ordinary geom face
259    *
260    * It's parametrization is within [0,1] range.
261    * It redefines Adaptor3d_Surface::Value(U,V) where U and V are within [0,1]
262    */
263   // --------------------------------------------------------------------
264   class TSideFace: public Adaptor3d_Surface
265   {
266     typedef boost::shared_ptr<BRepAdaptor_Surface> PSurface;
267
268     int                             myID; //!< in-block ID
269     // map used to find out real UV by it's normalized UV
270     TParam2ColumnMap*               myParamToColumnMap;
271     PSurface                        mySurface;
272     TopoDS_Edge                     myBaseEdge;
273     std::map< int, PSurface >       myShapeID2Surf;
274     // first and last normalized params and orientaion for each component or it-self
275     std::vector< std::pair< double, double> > myParams; // select my columns in myParamToColumnMap
276     bool                            myIsForward;
277     std::vector< TSideFace* >       myComponents;
278     SMESH_MesherHelper              myHelper;
279   public:
280     TSideFace( SMESH_Mesh&                mesh,
281                const int                  faceID,
282                const Prism_3D::TQuadList& quadList,
283                const TopoDS_Edge&         baseEdge,
284                TParam2ColumnMap*          columnsMap,
285                const double               first = 0.0,
286                const double               last  = 1.0);
287     TSideFace( SMESH_Mesh&                                       mesh,
288                const std::vector< TSideFace* >&                  components,
289                const std::vector< std::pair< double, double> > & params);
290     TSideFace( const TSideFace& other );
291     ~TSideFace();
292     bool IsComplex() const
293     { return ( NbComponents() > 0 || myParams[0].first != 0. || myParams[0].second != 1. ); }
294     int FaceID() const { return myID; }
295     SMESH_Mesh* GetMesh() const { return myHelper.GetMesh(); }
296     TParam2ColumnMap* GetColumns() const { return myParamToColumnMap; }
297     gp_XY GetNodeUV(const TopoDS_Face& F, const SMDS_MeshNode* n, const SMDS_MeshNode* n2=0) const
298     { return ((SMESH_MesherHelper&) myHelper).SetSubShape(F), myHelper.GetNodeUV( F, n, n2 ); }
299     const TopoDS_Edge & BaseEdge() const { return myBaseEdge; }
300     int ColumnHeight() const {
301       if ( NbComponents() ) return GetComponent(0)->GetColumns()->begin()->second.size();
302       else                  return GetColumns()->begin()->second.size(); }
303     double GetColumns(const double U, TParam2ColumnIt & col1, TParam2ColumnIt& col2 ) const;
304     void GetNodesAtZ(const int Z, std::map<double, const SMDS_MeshNode* >& nodes ) const;
305     int NbComponents() const { return myComponents.size(); }
306     TSideFace* GetComponent(const int i) const { return myComponents.at( i ); }
307     void SetComponent(const int i, TSideFace* c)
308     { if ( myComponents[i] ) delete myComponents[i]; myComponents[i]=c; }
309     TSideFace* GetComponent(const double U, double& localU) const;
310     bool IsForward() const { return myIsForward; }
311     // boundary geometry for a face
312     Adaptor3d_Surface* Surface() const { return new TSideFace( *this ); }
313     bool GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const;
314     Adaptor2d_Curve2d* HorizPCurve(const bool isTop, const TopoDS_Face& horFace) const;
315     Adaptor3d_Curve* HorizCurve(const bool isTop) const;
316     Adaptor3d_Curve* VertiCurve(const bool isMax) const;
317     TopoDS_Edge GetEdge( const int edge ) const;
318     int InsertSubShapes( TBlockShapes& shapeMap ) const;
319     // redefine Adaptor methods
320     gp_Pnt Value(const Standard_Real U,const Standard_Real V) const;
321     // debug
322     void dumpNodes(int nbNodes) const;
323   };
324
325   // --------------------------------------------------------------------
326   /*!
327    * \brief Class emulating geometry of a vertical edge
328    */
329   // --------------------------------------------------------------------
330   class STDMESHERS_EXPORT TVerticalEdgeAdaptor: public Adaptor3d_Curve
331   {
332     const TNodeColumn* myNodeColumn;
333   public:
334     TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter );
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 geometry of a hirizontal edge
345    */
346   // --------------------------------------------------------------------
347   class STDMESHERS_EXPORT THorizontalEdgeAdaptor: public Adaptor3d_Curve
348   {
349     const TSideFace* mySide;
350     double           myV;
351   public:
352     THorizontalEdgeAdaptor( const TSideFace* sideFace, const bool isTop)
353       :mySide(sideFace), myV( isTop ? 1.0 : 0.0 ) {}
354     gp_Pnt Value(const Standard_Real U) const;
355     Standard_Real FirstParameter() const { return 0; }
356     Standard_Real LastParameter() const { return 1; }
357     // debug
358     void dumpNodes(int nbNodes) const;
359   };
360
361   // --------------------------------------------------------------------
362   /*!
363    * \brief Class emulating pcurve on a hirizontal face
364    */
365   // --------------------------------------------------------------------
366   class STDMESHERS_EXPORT TPCurveOnHorFaceAdaptor: public Adaptor2d_Curve2d
367   {
368     std::map< double, gp_XY > myUVmap; // normalized parameter to UV on a horizontal face
369   public:
370     TPCurveOnHorFaceAdaptor( const TSideFace*   sideFace,
371                              const bool         isTop,
372                              const TopoDS_Face& horFace);
373     gp_Pnt2d Value(const Standard_Real U) const;
374     Standard_Real FirstParameter() const { return 0; }
375     Standard_Real LastParameter() const { return 1; }
376   };
377
378   bool                  myNotQuadOnTop;
379   SMESH_MesherHelper*   myHelper;
380   TBlockShapes          myShapeIDMap;
381   SMESH_ComputeErrorPtr myError;
382
383   // container of 4 side faces
384   TSideFace*            mySide; 
385   // node columns for each base edge
386   std::vector< TParam2ColumnMap >                       myParam2ColumnMaps;
387   // to find a column for a node by edge SMESHDS Index
388   std::map< int, std::pair< TParam2ColumnMap*, bool > > myShapeIndex2ColumnMap;
389
390   /*!
391    * \brief store error and comment and then return ( error == COMPERR_OK )
392    */
393   bool error(int error, const SMESH_Comment& comment = "") {
394     myError = SMESH_ComputeError::New(error,comment);
395     return myError->IsOK();
396   }
397   /*!
398    * \brief Prints a script creating a normal grid on the prism side
399    */
400   void faceGridToPythonDump(const SMESH_Block::TShapeID face,
401                             const int                   nb=10);
402
403 }; // class StdMeshers_PrismAsBlock
404
405 // ===============================================
406 /*!
407  * \brief Tool building internal nodes in a prism
408  */
409 struct StdMeshers_Sweeper
410 {
411   std::vector< TNodeColumn* > myBndColumns; // boundary nodes
412   std::vector< TNodeColumn* > myIntColumns; // internal nodes
413
414   bool ComputeNodes( SMESH_MesherHelper& helper,
415                      const double        tol,
416                      const bool          allowHighBndError );
417
418 private:
419
420   gp_XYZ bndPoint( int iP, int z ) const
421   { return SMESH_TNodeXYZ( (*myBndColumns[ iP ])[ z ]); }
422
423   gp_XYZ intPoint( int iP, int z ) const
424   { return SMESH_TNodeXYZ( (*myIntColumns[ iP ])[ z ]); }
425
426   static bool projectIntPoints(const std::vector< gp_XYZ >& fromBndPoints,
427                                const std::vector< gp_XYZ >& toBndPoints,
428                                const std::vector< gp_XYZ >& fromIntPoints,
429                                std::vector< gp_XYZ >&       toIntPoints,
430                                StdMeshers_ProjectionUtils::TrsfFinder3D& trsf,
431                                std::vector< gp_XYZ > *      bndError);
432
433   static void applyBoundaryError(const std::vector< gp_XYZ >& bndPoints,
434                                  const std::vector< gp_XYZ >& bndError1,
435                                  const std::vector< gp_XYZ >& bndError2,
436                                  const double                 r,
437                                  std::vector< gp_XYZ >&       toIntPoints,
438                                  std::vector< double >&       int2BndDist);
439 };
440
441 // ===============================================
442 /*!
443  * \brief Algo building prisms on a prism shape
444  */
445 class STDMESHERS_EXPORT StdMeshers_Prism_3D: public SMESH_3D_Algo
446 {
447 public:
448   StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen);
449   virtual ~StdMeshers_Prism_3D();
450
451   virtual bool CheckHypothesis(SMESH_Mesh&                          aMesh,
452                                const TopoDS_Shape&                  aShape,
453                                SMESH_Hypothesis::Hypothesis_Status& aStatus);
454
455   virtual bool Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape);
456
457   virtual bool Evaluate(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape,
458                         MapShapeNbElems& aResMap);
459
460   /*!
461    * \brief Enable removal of quadrangles from the bottom face and
462    * triangles creation there by projection from the top
463    * (sole face meshed with triangles is considered to be a bottom one).
464    * If there are two faces with triangles, triangles must
465    * be of the same topology, else the algo fails.
466    * The method must be called before Compute()
467    */
468   void ProjectTriangles() { myProjectTriangles = true; }
469
470   /*!
471    * \brief Create prisms
472     * \param nodeColumns - columns of nodes generated from nodes of a mesh face
473     * \param helper - helper initialized by mesh and shape to add prisms to
474    */
475   static bool AddPrisms( std::vector<const TNodeColumn*> & nodeColumns,
476                          SMESH_MesherHelper*               helper);
477
478   static bool IsApplicable(const TopoDS_Shape & aShape, bool toCheckAll);
479
480  private:
481
482   /*!
483    * \brief Analyse shape geometry and mesh.
484     * If there are triangles on one of faces, it becomes 'bottom'
485    */
486   bool initPrism(Prism_3D::TPrismTopo& thePrism,
487                  const TopoDS_Shape&   theSolid,
488                  const bool            selectBottom = true);
489
490   /*!
491    * \brief Fill thePrism.myWallQuads and thePrism.myTopEdges
492    */
493   bool getWallFaces( Prism_3D::TPrismTopo& thePrism,
494                      const int             totalNbFaces);
495
496   /*!
497    * \brief Compute mesh on a SOLID
498    */
499   bool compute(const Prism_3D::TPrismTopo& thePrism);
500
501   /*!
502    * \brief Compute 2D mesh on walls FACEs of a prism
503    */
504   bool computeWalls(const Prism_3D::TPrismTopo& thePrism);
505
506   /*!
507    * \brief Returns a source EDGE of propagation to a given EDGE
508    */
509   TopoDS_Edge findPropagationSource( const TopoDS_Edge& E );
510
511   /*!
512    * \brief Find correspondence between bottom and top nodes.
513    *  If elements on the bottom and top faces are topologically different,
514    *  and projection is possible and allowed, perform the projection
515     * \retval bool - is a success or not
516    */
517   bool assocOrProjBottom2Top( const gp_Trsf & bottomToTopTrsf,
518                               const Prism_3D::TPrismTopo& thePrism);
519
520   /*!
521    * \brief Remove quadrangles from the top face and
522    *        create triangles there by projection from the bottom
523     * \retval bool - a success or not
524    */
525   bool projectBottomToTop( const gp_Trsf & bottomToTopTrsf,
526                            const Prism_3D::TPrismTopo& thePrism );
527
528   /*!
529    * \brief Compute tolerance to pass to StdMeshers_Sweeper
530    */
531   double getSweepTolerance( const Prism_3D::TPrismTopo& thePrism );
532
533   /*!
534    * \brief Defines if it's safe to use the block approach
535    */
536   bool isSimpleBottom( const Prism_3D::TPrismTopo& thePrism );
537
538   /*!
539    * \brief Project mesh faces from a source FACE of one prism to
540    *        a source FACE of another prism
541    *  \retval bool - a success or not
542    */
543   bool project2dMesh(const TopoDS_Face& source, const TopoDS_Face& target);
544
545   /*!
546    * \brief Set projection coordinates of a node to a face and it's sub-shapes
547     * \param faceID - the face given by in-block ID
548     * \param params - node normalized parameters
549     * \retval bool - is a success
550    */
551   bool setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z );
552
553   /*!
554    * \brief If (!isOK), sets the error to a sub-mesh of a current SOLID
555    */
556   bool toSM( bool isOK );
557
558   /*!
559    * \brief Return index of a shape
560    */
561   int shapeID( const TopoDS_Shape& S );
562
563 private:
564
565   bool myProjectTriangles;
566   bool mySetErrorToSM;
567   bool myUseBlock;
568
569   StdMeshers_PrismAsBlock myBlock;
570   SMESH_MesherHelper*     myHelper;
571
572   std::vector<gp_XYZ>     myShapeXYZ; // point on each sub-shape of the block
573
574   // map of bottom nodes to the column of nodes above them
575   // (the column includes the bottom node)
576   TNode2ColumnMap         myBotToColumnMap;
577
578   TopTools_IndexedMapOfShape* myPropagChains;
579
580 }; // class StdMeshers_Prism_3D
581
582 #endif