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