SALOME - SMESH
StdMeshers_Prism_3D.hxx
Go to the documentation of this file.
1 // Copyright (C) 2007-2008 CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 // SMESH SMESH : implementaion of SMESH idl descriptions
23 // File : StdMeshers_Prism_3D.hxx
24 // Module : SMESH
25 //
26 #ifndef _SMESH_Prism_3D_HXX_
27 #define _SMESH_Prism_3D_HXX_
28 
29 #include "SMESH_StdMeshers.hxx"
30 
31 #include "SMESH_3D_Algo.hxx"
32 #include "SMDS_TypeOfPosition.hxx"
33 #include "SMDS_MeshNode.hxx"
34 #include "SMESH_Block.hxx"
35 #include "SMESH_Mesh.hxx"
36 #include "SMESHDS_Mesh.hxx"
37 #include "SMESH_subMesh.hxx"
38 #include "SMESH_MesherHelper.hxx"
39 #include "SMESH_Comment.hxx"
40 
41 #include <vector>
42 
43 #include <Adaptor3d_Curve.hxx>
44 #include <Adaptor3d_Surface.hxx>
45 #include <Adaptor2d_Curve2d.hxx>
46 #include <BRepAdaptor_Surface.hxx>
47 #include <TopTools_IndexedMapOfOrientedShape.hxx>
48 #include <gp_XYZ.hxx>
49 
50 
51 class SMESHDS_SubMesh;
52 class TopoDS_Edge;
53 class TopoDS_Faces;
54 struct TNode;
55 
56 //typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*> TNodeNodeMap;
57 typedef std::vector<const SMDS_MeshNode* > TNodeColumn;
58 
59 // map of bottom nodes to the column of nodes above them
60 // (the column includes the bottom nodes)
61 typedef std::map< TNode, TNodeColumn > TNode2ColumnMap;
62 typedef std::map< double, TNodeColumn > TParam2ColumnMap;
63 typedef std::map< double, TNodeColumn >::const_iterator TParam2ColumnIt;
64 
65 typedef TopTools_IndexedMapOfOrientedShape TBlockShapes;
66 
67 // ===============================================
71 // ===============================================
72 
73 struct TNode
74 {
76  gp_XYZ myParams;
77 
78  gp_XYZ GetCoords() const { return gp_XYZ( myNode->X(), myNode->Y(), myNode->Z() ); }
79  gp_XYZ GetParams() const { return myParams; }
80  gp_XYZ& ChangeParams() { return myParams; }
81  bool HasParams() const { return myParams.X() >= 0.0; }
83  { return myNode ? myNode->GetPosition()->GetTypeOfPosition() : SMDS_TOP_UNSPEC; }
84  bool IsNeighbor( const TNode& other ) const;
85 
86  TNode(const SMDS_MeshNode* node = 0): myNode(node), myParams(-1,-1,-1) {}
87  bool operator < (const TNode& other) const { return myNode->GetID() < other.myNode->GetID(); }
88 };
89 
90 // ===============================================================
97 // ===============================================================
98 
100 {
101 public:
106 
108 
118  bool Init(SMESH_MesherHelper* helper, const TopoDS_Shape& shape3D);
119 
123  SMESH_ComputeErrorPtr GetError() const { return myError; }
124 
129  int VerticalSize() const { return myParam2ColumnMaps[0].begin()->second.size(); }
130 
131  bool HasNotQuadElemOnTop() const { return myNotQuadOnTop; }
132 
138  const TNodeColumn* GetNodeColumn(const SMDS_MeshNode* node) const;
139 
146  const TParam2ColumnMap& GetParam2ColumnMap(const int baseEdgeID,
147  bool & isReverse)
148  {
149  std::pair< TParam2ColumnMap*, bool > & col_frw =
150  myShapeIndex2ColumnMap[ baseEdgeID ];
151  isReverse = !col_frw.second;
152  return * col_frw.first;
153  }
154 
159  SMESH_Mesh* Mesh() const { return myHelper->GetMesh(); }
160 
165  SMESHDS_Mesh* MeshDS() const { return Mesh()->GetMeshDS(); }
166 
172  SMESH_subMesh* SubMesh(const int shapeID) const
173  { return Mesh()->GetSubMesh( Shape( shapeID )); }
174 
180  SMESHDS_SubMesh* SubMeshDS(const int shapeID) const
181  { return SubMesh(shapeID)->GetSubMeshDS(); }
182 
188  const TopoDS_Shape& Shape(const int shapeID) const
189  { return myShapeIDMap( shapeID ); }
190 
196  int ShapeID(const TopoDS_Shape& shape) const
197  { return myShapeIDMap.FindIndex( shape ); }
198 
207  static bool IsForwardEdge(SMESHDS_Mesh* meshDS,
208  const TParam2ColumnMap& columnsMap,
209  const TopoDS_Edge & bottomEdge,
210  const int sideFaceID);
219  static bool GetWallFaces( SMESH_Mesh* mesh,
220  const TopoDS_Shape & mainShape,
221  const TopoDS_Shape & bottomFace,
222  const std::list< TopoDS_Edge >& bottomEdges,
223  std::list< TopoDS_Face >& wallFaces);
224 
225 private:
226 
227  // --------------------------------------------------------------------
235  // --------------------------------------------------------------------
237  {
238  int myID;
239  // map used to find out real UV by it's normalized UV
241  BRepAdaptor_Surface mySurface;
242  TopoDS_Edge myBaseEdge;
243  // first and last normalized params and orientaion for each component or it-self
244  std::vector< std::pair< double, double> > myParams;
246  std::vector< TSideFace* > myComponents;
248  public:
250  const int faceID,
251  const TopoDS_Face& face,
252  const TopoDS_Edge& baseEdge,
253  TParam2ColumnMap* columnsMap,
254  const double first = 0.0,
255  const double last = 1.0);
256  TSideFace( const std::vector< TSideFace* >& components,
257  const std::vector< std::pair< double, double> > & params);
258  TSideFace( const TSideFace& other );
260  bool IsComplex() const
261  { return ( NbComponents() > 0 || myParams[0].first != 0. || myParams[0].second != 1. ); }
262  int FaceID() const { return myID; }
263  TParam2ColumnMap* GetColumns() const { return myParamToColumnMap; }
264  gp_XY GetNodeUV(const TopoDS_Face& F, const SMDS_MeshNode* n) const
265  { return myHelper->GetNodeUV( F, n ); }
266  const TopoDS_Edge & BaseEdge() const { return myBaseEdge; }
267  int ColumnHeight() const {
268  if ( NbComponents() ) return GetComponent(0)->GetColumns()->begin()->second.size();
269  else return GetColumns()->begin()->second.size(); }
270  double GetColumns(const double U, TParam2ColumnIt & col1, TParam2ColumnIt& col2 ) const;
271  int NbComponents() const { return myComponents.size(); }
272  TSideFace* GetComponent(const int i) const { return myComponents.at( i ); }
273  void SetComponent(const int i, TSideFace* c)
274  { if ( myComponents[i] ) delete myComponents[i]; myComponents[i]=c; }
275  TSideFace* GetComponent(const double U, double& localU) const;
276  bool IsForward() const { return myIsForward; }
277  // boundary geometry for a face
278  Adaptor3d_Surface* Surface() const { return new TSideFace( *this ); }
279  bool GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const;
280  Adaptor2d_Curve2d* HorizPCurve(const bool isTop, const TopoDS_Face& horFace) const;
281  Adaptor3d_Curve* HorizCurve(const bool isTop) const;
282  Adaptor3d_Curve* VertiCurve(const bool isMax) const;
283  TopoDS_Edge GetEdge( const int edge ) const;
284  int InsertSubShapes( TBlockShapes& shapeMap ) const;
285  // redefine Adaptor methods
286  gp_Pnt Value(const Standard_Real U,const Standard_Real V) const;
287  };
288 
289  // --------------------------------------------------------------------
293  // --------------------------------------------------------------------
295  {
297  public:
298  TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter );
299  gp_Pnt Value(const Standard_Real U) const;
300  Standard_Real FirstParameter() const { return 0; }
301  Standard_Real LastParameter() const { return 1; }
302  };
303 
304  // --------------------------------------------------------------------
308  // --------------------------------------------------------------------
310  {
312  double myV;
313  public:
314  THorizontalEdgeAdaptor( const TSideFace* sideFace, const bool isTop)
315  :mySide(sideFace), myV( isTop ? 1.0 : 0.0 ) {}
316  gp_Pnt Value(const Standard_Real U) const;
317  Standard_Real FirstParameter() const { return 0; }
318  Standard_Real LastParameter() const { return 1; }
319  };
320 
321  // --------------------------------------------------------------------
325  // --------------------------------------------------------------------
327  {
329  int myZ;
330  TopoDS_Face myFace;
331  public:
333  const bool isTop,
334  const TopoDS_Face& horFace)
335  : mySide(sideFace), myFace(horFace), myZ(isTop ? mySide->ColumnHeight() - 1 : 0 ) {}
336  gp_Pnt2d Value(const Standard_Real U) const;
337  Standard_Real FirstParameter() const { return 0; }
338  Standard_Real LastParameter() const { return 1; }
339  };
340  // --------------------------------------------------------------------
341 
345 
346  // container of 4 side faces
348  // node columns for each base edge
349  std::vector< TParam2ColumnMap > myParam2ColumnMaps;
350  // to find a column for a node by edge SMESHDS Index
351  std::map< int, std::pair< TParam2ColumnMap*, bool > > myShapeIndex2ColumnMap;
352 
357  bool error(int error, const SMESH_Comment& comment = "") {
358  myError = SMESH_ComputeError::New(error,comment);
359  return myError->IsOK();
360  }
361  //std::vector< SMESH_subMesh* > mySubMeshesVec; // submesh by in-block id
362 };
363 
364 // =============================================
368 // =============================================
369 
371 {
372 public:
373  StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen);
375 
376  virtual bool CheckHypothesis(SMESH_Mesh& aMesh,
377  const TopoDS_Shape& aShape,
379 
380  virtual bool Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape);
381 
382  virtual bool Evaluate(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape,
383  MapShapeNbElems& aResMap);
384 
393  void ProjectTriangles() { myProjectTriangles = true; }
394 
400  static void AddPrisms( std::vector<const TNodeColumn*> & nodeColumns,
401  SMESH_MesherHelper* helper);
402 
403 private:
404 
412 
419 
426  bool setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z );
427 
428 private:
429 
431 
434 
435  std::vector<gp_XYZ> myShapeXYZ; // point on each sub-shape
436 
437  // map of bottom nodes to the column of nodes above them
438  // (the column includes the bottom node)
440 };
441 
442 #endif
SMDS_TypeOfPosition
@ SMDS_TOP_UNSPEC
std::map< SMESH_subMesh *, std::vector< int > > MapShapeNbElems
Definition: SMESH_Algo.hxx:55
boost::shared_ptr< SMESH_ComputeError > SMESH_ComputeErrorPtr
std::map< double, TNodeColumn > TParam2ColumnMap
std::vector< const SMDS_MeshNode * > TNodeColumn
It helps meshers to add elements.
#define STDMESHERS_EXPORT
std::map< double, TNodeColumn > TParam2ColumnMap
TopTools_IndexedMapOfOrientedShape TBlockShapes
std::map< TNode, TNodeColumn > TNode2ColumnMap
std::map< double, TNodeColumn >::const_iterator TParam2ColumnIt
std::vector< const SMDS_MeshNode * > TNodeColumn
int GetID() const
double X() const
double Z() const
double Y() const
const SMDS_PositionPtr & GetPosition() const
Class to generate string from any type.
gp_XY GetNodeUV(const TopoDS_Face &F, const SMDS_MeshNode *n, const SMDS_MeshNode *inFaceNode=0, bool *check=0) const
Return node UV on face.
Class emulating geometry of a hirizontal edge.
THorizontalEdgeAdaptor(const TSideFace *sideFace, const bool isTop)
gp_Pnt Value(const Standard_Real U) const
Class emulating pcurve on a hirizontal face.
gp_Pnt2d Value(const Standard_Real U) const
TPCurveOnHorFaceAdaptor(const TSideFace *sideFace, const bool isTop, const TopoDS_Face &horFace)
Class representing a part of a geom face or a union of seleral faces. Or just an ordinary geom face.
const TopoDS_Edge & BaseEdge() const
TSideFace * GetComponent(const double U, double &localU) const
double GetColumns(const double U, TParam2ColumnIt &col1, TParam2ColumnIt &col2) const
bool GetPCurves(Adaptor2d_Curve2d *pcurv[4]) const
TSideFace * GetComponent(const int i) const
gp_XY GetNodeUV(const TopoDS_Face &F, const SMDS_MeshNode *n) const
Adaptor3d_Curve * VertiCurve(const bool isMax) const
Adaptor3d_Surface * Surface() const
TParam2ColumnMap * GetColumns() const
Adaptor2d_Curve2d * HorizPCurve(const bool isTop, const TopoDS_Face &horFace) const
TopoDS_Edge GetEdge(const int edge) const
std::vector< std::pair< double, double > > myParams
Adaptor3d_Curve * HorizCurve(const bool isTop) const
TSideFace(const TSideFace &other)
void SetComponent(const int i, TSideFace *c)
TSideFace(SMESH_MesherHelper *helper, const int faceID, const TopoDS_Face &face, const TopoDS_Edge &baseEdge, TParam2ColumnMap *columnsMap, const double first=0.0, const double last=1.0)
int InsertSubShapes(TBlockShapes &shapeMap) const
TSideFace(const std::vector< TSideFace * > &components, const std::vector< std::pair< double, double > > &params)
gp_Pnt Value(const Standard_Real U, const Standard_Real V) const
std::vector< TSideFace * > myComponents
Class emulating geometry of a vertical edge.
gp_Pnt Value(const Standard_Real U) const
TVerticalEdgeAdaptor(const TParam2ColumnMap *columnsMap, const double parameter)
Tool analyzing and giving access to a prism geometry treating it like a block, i.e....
const TNodeColumn * GetNodeColumn(const SMDS_MeshNode *node) const
Return pointer to column of nodes.
SMESH_ComputeErrorPtr GetError() const
Return problem description.
std::map< int, std::pair< TParam2ColumnMap *, bool > > myShapeIndex2ColumnMap
static bool IsForwardEdge(SMESHDS_Mesh *meshDS, const TParam2ColumnMap &columnsMap, const TopoDS_Edge &bottomEdge, const int sideFaceID)
Check curve orientation of a bootom edge.
StdMeshers_PrismAsBlock()
Constructor. Initialization is needed.
SMESHDS_SubMesh * SubMeshDS(const int shapeID) const
Return submesh DS of a shape.
SMESH_Mesh * Mesh() const
Return pointer to mesh.
std::vector< TParam2ColumnMap > myParam2ColumnMaps
SMESH_subMesh * SubMesh(const int shapeID) const
Return submesh of a shape.
bool error(int error, const SMESH_Comment &comment="")
store error and comment and then return ( error == COMPERR_OK )
int ShapeID(const TopoDS_Shape &shape) const
Return in-block ID of a shape.
SMESH_MesherHelper * myHelper
const TParam2ColumnMap & GetParam2ColumnMap(const int baseEdgeID, bool &isReverse)
Return TParam2ColumnMap for a base edge.
SMESHDS_Mesh * MeshDS() const
Return pointer to mesh DS.
const TopoDS_Shape & Shape(const int shapeID) const
Return a in-block shape.
bool Init(SMESH_MesherHelper *helper, const TopoDS_Shape &shape3D)
Initialization.
SMESH_ComputeErrorPtr myError
static bool GetWallFaces(SMESH_Mesh *mesh, const TopoDS_Shape &mainShape, const TopoDS_Shape &bottomFace, const std::list< TopoDS_Edge > &bottomEdges, std::list< TopoDS_Face > &wallFaces)
Find wall faces by bottom edges.
int VerticalSize() const
Return number of nodes on every vertical edge.
Algo building prisms on a prism shape.
bool setFaceAndEdgesXYZ(const int faceID, const gp_XYZ &params, int z)
Set projection coordinates of a node to a face and it's subshapes.
StdMeshers_PrismAsBlock myBlock
SMESH_MesherHelper * myHelper
TNode2ColumnMap myBotToColumnMap
virtual bool Compute(SMESH_Mesh &aMesh, const TopoDS_Shape &aShape)
Computes mesh on a shape.
void ProjectTriangles()
Enable removal of quadrangles from the bottom face and triangles creation there by projection from th...
static void AddPrisms(std::vector< const TNodeColumn * > &nodeColumns, SMESH_MesherHelper *helper)
Create prisms.
StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen *gen)
virtual bool Evaluate(SMESH_Mesh &aMesh, const TopoDS_Shape &aShape, MapShapeNbElems &aResMap)
evaluates size of prospective mesh on a shape
virtual ~StdMeshers_Prism_3D()
virtual bool CheckHypothesis(SMESH_Mesh &aMesh, const TopoDS_Shape &aShape, SMESH_Hypothesis::Hypothesis_Status &aStatus)
Check hypothesis definition to mesh a shape.
std::vector< gp_XYZ > myShapeXYZ
bool assocOrProjBottom2Top()
Find correspondence between bottom and top nodes. If elements on the bottom and top faces are topolog...
bool projectBottomToTop()
Remove quadrangles from the top face and create triangles there by projection from the bottom.
static SMESH_ComputeErrorPtr New(int error=COMPERR_OK, std::string comment="", const SMESH_Algo *algo=0)
Structure containing node relative data.
gp_XYZ GetParams() const
bool IsNeighbor(const TNode &other) const
gp_XYZ GetCoords() const
bool HasParams() const
const SMDS_MeshNode * myNode
bool operator<(const TNode &other) const
gp_XYZ & ChangeParams()
SMDS_TypeOfPosition GetPositionType() const
TNode(const SMDS_MeshNode *node=0)