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