Version: 6.3.1

src/StdMeshers/StdMeshers_Prism_3D.hxx

Go to the documentation of this file.
00001 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
00002 //
00003 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
00004 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
00005 //
00006 // This library is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 2.1 of the License.
00010 //
00011 // This library is distributed in the hope that it will be useful,
00012 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00013 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00014 // Lesser General Public License for more details.
00015 //
00016 // You should have received a copy of the GNU Lesser General Public
00017 // License along with this library; if not, write to the Free Software
00018 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00019 //
00020 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00021 //
00022 
00023 //  SMESH SMESH : implementaion of SMESH idl descriptions
00024 //  File   : StdMeshers_Prism_3D.hxx
00025 //  Module : SMESH
00026 //
00027 #ifndef _SMESH_Prism_3D_HXX_
00028 #define _SMESH_Prism_3D_HXX_
00029 
00030 #include "SMESH_StdMeshers.hxx"
00031 
00032 #include "SMESH_3D_Algo.hxx"
00033 #include "SMDS_TypeOfPosition.hxx"
00034 #include "SMDS_MeshNode.hxx"
00035 #include "SMESH_Block.hxx"
00036 #include "SMESH_Mesh.hxx"
00037 #include "SMESHDS_Mesh.hxx"
00038 #include "SMESH_subMesh.hxx"
00039 #include "SMESH_MesherHelper.hxx"
00040 #include "SMESH_Comment.hxx"
00041 
00042 #include <vector>
00043 
00044 #include <Adaptor3d_Curve.hxx>
00045 #include <Adaptor3d_Surface.hxx>
00046 #include <Adaptor2d_Curve2d.hxx>
00047 #include <BRepAdaptor_Surface.hxx>
00048 #include <TopTools_IndexedMapOfOrientedShape.hxx>
00049 #include <gp_XYZ.hxx>
00050 #include <gp_Trsf.hxx>
00051 
00052 
00053 class SMESHDS_SubMesh;
00054 class TopoDS_Edge;
00055 class TopoDS_Faces;
00056 struct TNode;
00057 
00058 //typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*> TNodeNodeMap;
00059 typedef std::vector<const SMDS_MeshNode* > TNodeColumn;
00060 
00061 // map of bottom nodes to the column of nodes above them
00062 // (the column includes the bottom nodes)
00063 typedef std::map< TNode, TNodeColumn >  TNode2ColumnMap;
00064 typedef std::map< double, TNodeColumn > TParam2ColumnMap;
00065 typedef std::map< double, TNodeColumn >::const_iterator TParam2ColumnIt;
00066 
00067 typedef TopTools_IndexedMapOfOrientedShape TBlockShapes;
00068 
00069 // ===============================================
00073 // ===============================================
00074 
00075 struct TNode
00076 {
00077   const SMDS_MeshNode* myNode;
00078   mutable gp_XYZ       myParams;
00079 
00080   gp_XYZ GetCoords() const { return gp_XYZ( myNode->X(), myNode->Y(), myNode->Z() ); }
00081   gp_XYZ GetParams() const { return myParams; }
00082   gp_XYZ& ChangeParams() const { return myParams; }
00083   bool HasParams() const { return myParams.X() >= 0.0; }
00084   SMDS_TypeOfPosition GetPositionType() const
00085   { return myNode ? myNode->GetPosition()->GetTypeOfPosition() : SMDS_TOP_UNSPEC; }
00086   bool IsNeighbor( const TNode& other ) const;
00087 
00088   TNode(const SMDS_MeshNode* node = 0): myNode(node), myParams(-1,-1,-1) {}
00089   bool operator < (const TNode& other) const { return myNode->GetID() < other.myNode->GetID(); }
00090 };
00091 
00092 // ===============================================================
00099 // ===============================================================
00100 
00101 class STDMESHERS_EXPORT StdMeshers_PrismAsBlock: public SMESH_Block
00102 {
00103 public:
00107   StdMeshers_PrismAsBlock();
00108 
00109   ~StdMeshers_PrismAsBlock();
00110 
00120   bool Init(SMESH_MesherHelper* helper, const TopoDS_Shape& shape3D);
00121 
00125   SMESH_ComputeErrorPtr GetError() const { return myError; }
00126 
00130   void Clear();
00131 
00136   int VerticalSize() const { return myParam2ColumnMaps[0].begin()->second.size(); }
00137 
00138   bool HasNotQuadElemOnTop() const { return myNotQuadOnTop; }
00139 
00145   const TNodeColumn* GetNodeColumn(const SMDS_MeshNode* node) const;
00146 
00153   const TParam2ColumnMap& GetParam2ColumnMap(const int baseEdgeID,
00154                                              bool &    isReverse) const
00155   {
00156     std::pair< TParam2ColumnMap*, bool > col_frw =
00157       myShapeIndex2ColumnMap.find( baseEdgeID )->second;
00158     isReverse = !col_frw.second;
00159     return * col_frw.first;
00160   }
00161 
00167   bool GetLayersTransformation(std::vector<gp_Trsf> & trsf) const;
00168   
00173   SMESH_Mesh* Mesh() const { return myHelper->GetMesh(); }
00174 
00179   SMESHDS_Mesh* MeshDS() const { return Mesh()->GetMeshDS(); }
00180 
00186   SMESH_subMesh* SubMesh(const int shapeID) const
00187   { return Mesh()->GetSubMesh( Shape( shapeID )); }
00188 
00194   SMESHDS_SubMesh* SubMeshDS(const int shapeID) const
00195   { return SubMesh(shapeID)->GetSubMeshDS(); }
00196 
00202   const TopoDS_Shape& Shape(const int shapeID) const
00203   { return myShapeIDMap( shapeID ); }
00204 
00210   int ShapeID(const TopoDS_Shape& shape) const
00211   { return myShapeIDMap.FindIndex( shape ); }
00212 
00221   static bool IsForwardEdge(SMESHDS_Mesh*           meshDS,
00222                             const TParam2ColumnMap& columnsMap,
00223                             const TopoDS_Edge &     bottomEdge,
00224                             const int               sideFaceID);
00233   bool GetWallFaces( SMESH_Mesh*               mesh,
00234                      const TopoDS_Shape &      mainShape,
00235                      const TopoDS_Shape &      bottomFace,
00236                      std::list< TopoDS_Edge >& bottomEdges,
00237                      std::list< int > &        nbEInW,
00238                      std::list< TopoDS_Face >& wallFaces);
00239 
00240 private:
00241 
00242   // --------------------------------------------------------------------
00250   // --------------------------------------------------------------------
00251   class TSideFace: public Adaptor3d_Surface
00252   {
00253     int                             myID; 
00254     // map used to find out real UV by it's normalized UV
00255     TParam2ColumnMap*               myParamToColumnMap;
00256     BRepAdaptor_Surface             mySurface;
00257     TopoDS_Edge                     myBaseEdge;
00258     // first and last normalized params and orientaion for each component or it-self
00259     std::vector< std::pair< double, double> > myParams;
00260     bool                            myIsForward;
00261     std::vector< TSideFace* >       myComponents;
00262     SMESH_MesherHelper *            myHelper;
00263   public:
00264     TSideFace( SMESH_MesherHelper* helper,
00265                const int           faceID,
00266                const TopoDS_Face&  face,
00267                const TopoDS_Edge&  baseEdge,
00268                TParam2ColumnMap*   columnsMap,
00269                const double        first = 0.0,
00270                const double        last = 1.0);
00271     TSideFace( const std::vector< TSideFace* >&             components,
00272                const std::vector< std::pair< double, double> > & params);
00273     TSideFace( const TSideFace& other );
00274     ~TSideFace();
00275     bool IsComplex() const
00276     { return ( NbComponents() > 0 || myParams[0].first != 0. || myParams[0].second != 1. ); }
00277     int FaceID() const { return myID; }
00278     TParam2ColumnMap* GetColumns() const { return myParamToColumnMap; }
00279     gp_XY GetNodeUV(const TopoDS_Face& F, const SMDS_MeshNode* n) const
00280     { return myHelper->GetNodeUV( F, n ); }
00281     const TopoDS_Edge & BaseEdge() const { return myBaseEdge; }
00282     int ColumnHeight() const {
00283       if ( NbComponents() ) return GetComponent(0)->GetColumns()->begin()->second.size();
00284       else                  return GetColumns()->begin()->second.size(); }
00285     double GetColumns(const double U, TParam2ColumnIt & col1, TParam2ColumnIt& col2 ) const;
00286     int NbComponents() const { return myComponents.size(); }
00287     TSideFace* GetComponent(const int i) const { return myComponents.at( i ); }
00288     void SetComponent(const int i, TSideFace* c)
00289     { if ( myComponents[i] ) delete myComponents[i]; myComponents[i]=c; }
00290     TSideFace* GetComponent(const double U, double& localU) const;
00291     bool IsForward() const { return myIsForward; }
00292     // boundary geometry for a face
00293     Adaptor3d_Surface* Surface() const { return new TSideFace( *this ); }
00294     bool GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const;
00295     Adaptor2d_Curve2d* HorizPCurve(const bool isTop, const TopoDS_Face& horFace) const;
00296     Adaptor3d_Curve* HorizCurve(const bool isTop) const;
00297     Adaptor3d_Curve* VertiCurve(const bool isMax) const;
00298     TopoDS_Edge GetEdge( const int edge ) const;
00299     int InsertSubShapes( TBlockShapes& shapeMap ) const;
00300     // redefine Adaptor methods
00301     gp_Pnt Value(const Standard_Real U,const Standard_Real V) const;
00302     // debug
00303     void dumpNodes(int nbNodes) const;
00304   };
00305 
00306   // --------------------------------------------------------------------
00310   // --------------------------------------------------------------------
00311   class STDMESHERS_EXPORT TVerticalEdgeAdaptor: public Adaptor3d_Curve
00312   {
00313     const TNodeColumn* myNodeColumn;
00314   public:
00315     TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter );
00316     gp_Pnt Value(const Standard_Real U) const;
00317     Standard_Real FirstParameter() const { return 0; }
00318     Standard_Real LastParameter() const { return 1; }
00319     // debug
00320     void dumpNodes(int nbNodes) const;
00321   };
00322 
00323   // --------------------------------------------------------------------
00327   // --------------------------------------------------------------------
00328   class STDMESHERS_EXPORT THorizontalEdgeAdaptor: public Adaptor3d_Curve
00329   {
00330     const TSideFace* mySide;
00331     double           myV;
00332   public:
00333     THorizontalEdgeAdaptor( const TSideFace* sideFace, const bool isTop)
00334       :mySide(sideFace), myV( isTop ? 1.0 : 0.0 ) {}
00335     gp_Pnt Value(const Standard_Real U) const;
00336     Standard_Real FirstParameter() const { return 0; }
00337     Standard_Real LastParameter() const { return 1; }
00338     // debug
00339     void dumpNodes(int nbNodes) const;
00340   };
00341 
00342   // --------------------------------------------------------------------
00346   // --------------------------------------------------------------------
00347   class STDMESHERS_EXPORT TPCurveOnHorFaceAdaptor: public Adaptor2d_Curve2d
00348   {
00349     const TSideFace*  mySide;
00350     int               myZ;
00351     TopoDS_Face       myFace;
00352   public:
00353     TPCurveOnHorFaceAdaptor( const TSideFace*   sideFace,
00354                              const bool         isTop,
00355                              const TopoDS_Face& horFace)
00356       : mySide(sideFace), myFace(horFace), myZ(isTop ? mySide->ColumnHeight() - 1 : 0 ) {}
00357     gp_Pnt2d Value(const Standard_Real U) const;
00358     Standard_Real FirstParameter() const { return 0; }
00359     Standard_Real LastParameter() const { return 1; }
00360   };
00361 
00362   bool                  myNotQuadOnTop;
00363   SMESH_MesherHelper*   myHelper;
00364   TBlockShapes          myShapeIDMap;
00365   SMESH_ComputeErrorPtr myError;
00366 
00367   // container of 4 side faces
00368   TSideFace*            mySide; 
00369   // node columns for each base edge
00370   std::vector< TParam2ColumnMap >                       myParam2ColumnMaps;
00371   // to find a column for a node by edge SMESHDS Index
00372   std::map< int, std::pair< TParam2ColumnMap*, bool > > myShapeIndex2ColumnMap;
00373 
00377   bool error(int error, const SMESH_Comment& comment = "") {
00378     myError = SMESH_ComputeError::New(error,comment);
00379     return myError->IsOK();
00380   }
00381 };
00382 
00383 // =============================================
00387 // =============================================
00388 
00389 class STDMESHERS_EXPORT StdMeshers_Prism_3D: public SMESH_3D_Algo
00390 {
00391 public:
00392   StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen);
00393   virtual ~StdMeshers_Prism_3D();
00394 
00395   virtual bool CheckHypothesis(SMESH_Mesh&                          aMesh,
00396                                const TopoDS_Shape&                  aShape,
00397                                SMESH_Hypothesis::Hypothesis_Status& aStatus);
00398 
00399   virtual bool Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape);
00400 
00401   virtual bool Evaluate(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape,
00402                         MapShapeNbElems& aResMap);
00403 
00412   void ProjectTriangles() { myProjectTriangles = true; }
00413 
00419   static void AddPrisms( std::vector<const TNodeColumn*> & nodeColumns,
00420                          SMESH_MesherHelper*          helper);
00421 
00422 private:
00423 
00430   bool assocOrProjBottom2Top();
00431 
00437   bool projectBottomToTop();
00438 
00445   bool setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z );
00446 
00447 private:
00448 
00449   bool myProjectTriangles;
00450 
00451   StdMeshers_PrismAsBlock myBlock;
00452   SMESH_MesherHelper*     myHelper;
00453 
00454   std::vector<gp_XYZ>     myShapeXYZ; // point on each sub-shape of the block
00455 
00456   // map of bottom nodes to the column of nodes above them
00457   // (the column includes the bottom node)
00458   TNode2ColumnMap         myBotToColumnMap;
00459 };
00460 
00461 #endif
Copyright © 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
Copyright © 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS