Version: 6.3.1

src/SMESH/SMESH_MesherHelper.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 // File:      SMESH_MesherHelper.hxx
00024 // Created:   15.02.06 14:48:09
00025 // Author:    Sergey KUUL
00026 //
00027 #ifndef SMESH_MesherHelper_HeaderFile
00028 #define SMESH_MesherHelper_HeaderFile
00029 
00030 #include "SMESH_SMESH.hxx"
00031 
00032 #include "SMESH_MeshEditor.hxx" // needed for many meshers
00033 #include <SMDS_MeshNode.hxx>
00034 #include <SMDS_QuadraticEdge.hxx>
00035 
00036 #include <Geom_Surface.hxx>
00037 #include <TopoDS_Face.hxx>
00038 #include <TopoDS_Shape.hxx>
00039 #include <gp_Pnt2d.hxx>
00040 
00041 #include <map>
00042 #include <vector>
00043 
00044 class GeomAPI_ProjectPointOnSurf;
00045 class GeomAPI_ProjectPointOnCurve;
00046 class SMESH_ProxyMesh;
00047 
00048 typedef std::map<SMESH_TLink, const SMDS_MeshNode*>           TLinkNodeMap;
00049 typedef std::map<SMESH_TLink, const SMDS_MeshNode*>::iterator ItTLinkNode;
00050 
00051 typedef SMDS_Iterator<const TopoDS_Shape*>  PShapeIterator;
00052 typedef boost::shared_ptr< PShapeIterator > PShapeIteratorPtr;
00053   
00054 typedef std::vector<const SMDS_MeshNode* > TNodeColumn;
00055 typedef std::map< double, TNodeColumn >    TParam2ColumnMap;
00056 
00057 typedef gp_XY (*xyFunPtr)(const gp_XY& uv1, const gp_XY& uv2);
00058 
00059 //=======================================================================
00068 //=======================================================================
00069 
00070 class SMESH_EXPORT SMESH_MesherHelper
00071 {
00072 public:
00073   // ---------- PUBLIC UTILITIES ----------
00074   
00081   static bool IsMedium(const SMDS_MeshNode*      node,
00082                        const SMDSAbs_ElementType typeToCheck = SMDSAbs_All);
00083 
00097   static bool LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
00098                               const TopoDS_Face& theFace,
00099                               const TopoDS_Edge& theBaseEdge,
00100                               SMESHDS_Mesh*      theMesh,
00101                               SMESH_ProxyMesh*   theProxyMesh=0);
00108   static TopoDS_Shape GetSubShapeByNode(const SMDS_MeshNode* node,
00109                                         const SMESHDS_Mesh*  meshDS);
00110 
00117   static int WrapIndex(const int ind, const int nbNodes) {
00118     if ( ind < 0 ) return nbNodes + ind % nbNodes;
00119     if ( ind >= nbNodes ) return ind % nbNodes;
00120     return ind;
00121   }
00122 
00126   static int NbAncestors(const TopoDS_Shape& shape,
00127                          const SMESH_Mesh&   mesh,
00128                          TopAbs_ShapeEnum    ancestorType=TopAbs_SHAPE);
00132   static PShapeIteratorPtr GetAncestors(const TopoDS_Shape& shape,
00133                                         const SMESH_Mesh&   mesh,
00134                                         TopAbs_ShapeEnum    ancestorType);
00135 
00139   static TopAbs_Orientation GetSubShapeOri(const TopoDS_Shape& shape,
00140                                            const TopoDS_Shape& subShape);
00141 
00142   static bool IsSubShape( const TopoDS_Shape& shape, const TopoDS_Shape& mainShape );
00143 
00144   static bool IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh );
00145 
00146   static double MaxTolerance( const TopoDS_Shape& shape );
00147 
00148   static bool IsClosedEdge( const TopoDS_Edge& anEdge );
00149 
00150   static TopoDS_Vertex IthVertex( const bool is2nd, TopoDS_Edge anEdge, const bool CumOri=true );
00151 
00152 
00153 public:
00154   // ---------- PUBLIC INSTANCE METHODS ----------
00155 
00156   // constructor
00157   SMESH_MesherHelper(SMESH_Mesh& theMesh);
00158 
00159   SMESH_Mesh* GetMesh() const { return myMesh; }
00160     
00161   SMESHDS_Mesh* GetMeshDS() const { return GetMesh()->GetMeshDS(); }
00162     
00167   bool IsQuadraticSubMesh(const TopoDS_Shape& theShape);
00171   void SetIsQuadratic(const bool theBuildQuadratic)
00172   { myCreateQuadratic = theBuildQuadratic; }
00176   bool GetIsQuadratic() const { return myCreateQuadratic; }
00177 
00182   void FixQuadraticElements(bool volumeOnly=true);
00183 
00189   void SetElementsOnShape(bool toSet) { mySetElemOnShape = toSet; }
00190 
00194   void SetSubShape(const int           subShapeID);
00195   void SetSubShape(const TopoDS_Shape& subShape);
00200   int GetSubShapeID() const { return myShapeID; }
00204   const TopoDS_Shape& GetSubShape() const  { return myShape; }
00205 
00209   SMDS_MeshNode* AddNode(double x, double y, double z, int ID = 0);
00213   SMDS_MeshEdge* AddEdge(const SMDS_MeshNode* n1,
00214                          const SMDS_MeshNode* n2,
00215                          const int id = 0, 
00216                          const bool force3d = true);
00220   SMDS_MeshFace* AddFace(const SMDS_MeshNode* n1,
00221                          const SMDS_MeshNode* n2,
00222                          const SMDS_MeshNode* n3,
00223                          const int id=0, 
00224                          const bool force3d = false);
00228   SMDS_MeshFace* AddFace(const SMDS_MeshNode* n1,
00229                          const SMDS_MeshNode* n2,
00230                          const SMDS_MeshNode* n3,
00231                          const SMDS_MeshNode* n4,
00232                          const int id = 0,
00233                          const bool force3d = false);
00234 
00238   SMDS_MeshFace* AddPolygonalFace (const std::vector<const SMDS_MeshNode*>& nodes,
00239                                    const int id = 0,
00240                                    const bool force3d = false);
00244   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00245                              const SMDS_MeshNode* n2,
00246                              const SMDS_MeshNode* n3,
00247                              const SMDS_MeshNode* n4,
00248                              const int id = 0,
00249                              const bool force3d = true);
00253   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00254                              const SMDS_MeshNode* n2,
00255                              const SMDS_MeshNode* n3,
00256                              const SMDS_MeshNode* n4,
00257                              const SMDS_MeshNode* n5,
00258                              const int id = 0,
00259                              const bool force3d = true);
00263   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00264                              const SMDS_MeshNode* n2,
00265                              const SMDS_MeshNode* n3,
00266                              const SMDS_MeshNode* n4,
00267                              const SMDS_MeshNode* n5,
00268                              const SMDS_MeshNode* n6,
00269                              const int id = 0, 
00270                              const bool force3d = true);
00274   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00275                              const SMDS_MeshNode* n2,
00276                              const SMDS_MeshNode* n3,
00277                              const SMDS_MeshNode* n4,
00278                              const SMDS_MeshNode* n5,
00279                              const SMDS_MeshNode* n6,
00280                              const SMDS_MeshNode* n7,
00281                              const SMDS_MeshNode* n8,
00282                              const int id = 0, 
00283                              bool force3d = true);
00284 
00288   SMDS_MeshVolume* AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
00289                                         const std::vector<int>&                  quantities,
00290                                         const int                                ID=0,
00291                                         const bool                               force3d = true);
00295   double GetNodeU(const TopoDS_Edge&   theEdge,
00296                   const SMDS_MeshNode* theNode,
00297                   const SMDS_MeshNode* inEdgeNode=0,
00298                   bool*                check=0);
00303   gp_XY GetNodeUV(const TopoDS_Face&   F,
00304                   const SMDS_MeshNode* n,
00305                   const SMDS_MeshNode* inFaceNode=0,
00306                   bool*                check=0) const;
00313   bool CheckNodeUV(const TopoDS_Face&   F,
00314                    const SMDS_MeshNode* n,
00315                    gp_XY&               uv,
00316                    const double         tol,
00317                    const bool           force=false,
00318                    double               distXYZ[4]=0) const;
00325   bool CheckNodeU(const TopoDS_Edge&   E,
00326                   const SMDS_MeshNode* n,
00327                   double&              u,
00328                   const double         tol,
00329                   const bool           force=false,
00330                   double               distXYZ[4]=0) const;
00334   static gp_XY GetMiddleUV(const Handle(Geom_Surface)& surface,
00335                            const gp_XY&                uv1,
00336                            const gp_XY&                uv2);
00344 #define gp_XY_FunPtr(meth) \
00345   static gp_XY __gpXY_##meth (const gp_XY& uv1, const gp_XY& uv2) { return uv1.meth( uv2 ); } \
00346   static xyFunPtr gp_XY_##meth = & __gpXY_##meth
00347 
00353   static gp_XY applyIn2D(const Handle(Geom_Surface)& surface,
00354                          const gp_XY&                uv1,
00355                          const gp_XY&                uv2,
00356                          xyFunPtr                    fun,
00357                          const bool                  resultInPeriod=true);
00358                           
00366   bool GetNodeUVneedInFaceNode(const TopoDS_Face& F = TopoDS_Face()) const;
00367 
00371   GeomAPI_ProjectPointOnSurf& GetProjector(const TopoDS_Face& F,
00372                                            TopLoc_Location&   loc,
00373                                            double             tol=0 ) const; 
00374 
00382   bool IsDegenShape(const int subShape) const
00383   { return myDegenShapeIds.find( subShape ) != myDegenShapeIds.end(); }
00389   bool HasDegeneratedEdges() const { return !myDegenShapeIds.empty(); }
00390 
00399   bool IsSeamShape(const int subShape) const
00400   { return mySeamShapeIds.find( subShape ) != mySeamShapeIds.end(); }
00409   bool IsSeamShape(const TopoDS_Shape& subShape) const
00410   { return IsSeamShape( GetMeshDS()->ShapeToIndex( subShape )); }
00415   bool IsRealSeam(const int subShape) const
00416   { return mySeamShapeIds.find( -subShape ) != mySeamShapeIds.end(); }
00421   bool IsRealSeam(const TopoDS_Shape& subShape) const
00422   { return IsRealSeam( GetMeshDS()->ShapeToIndex( subShape)); }
00428   bool HasSeam() const { return !mySeamShapeIds.empty(); }
00433   int GetPeriodicIndex() const { return myParIndex; }
00437   double GetOtherParam(const double param) const;
00438 
00445   const SMDS_MeshNode* GetMediumNode(const SMDS_MeshNode* n1,
00446                                      const SMDS_MeshNode* n2,
00447                                      const bool force3d);
00451   void AddTLinkNode(const SMDS_MeshNode* n1,
00452                     const SMDS_MeshNode* n2,
00453                     const SMDS_MeshNode* n12);
00457   void AddTLinkNodeMap(const TLinkNodeMap& aMap)
00458     { myTLinkNodeMap.insert(aMap.begin(), aMap.end()); }
00459 
00460   void AddTLinks(const SMDS_MeshEdge*   edge);
00461   void AddTLinks(const SMDS_MeshFace*   face);
00462   void AddTLinks(const SMDS_MeshVolume* vol);
00463 
00467   const TLinkNodeMap& GetTLinkNodeMap() const { return myTLinkNodeMap; }
00468 
00474   enum MType{ LINEAR, QUADRATIC, COMP };
00475   MType IsQuadraticMesh();
00476   
00477   virtual ~SMESH_MesherHelper();
00478 
00479 protected:
00480 
00487   gp_Pnt2d GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const;
00488 
00489   const SMDS_MeshNode* getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
00490                                                    const SMDS_MeshNode* n2,
00491                                                    bool                 force3d);
00492  private:
00493 
00494   // Forbiden copy constructor
00495   SMESH_MesherHelper (const SMESH_MesherHelper& theOther) {};
00496 
00497   // special map for using during creation of quadratic elements
00498   TLinkNodeMap    myTLinkNodeMap;
00499 
00500   std::set< int > myDegenShapeIds;
00501   std::set< int > mySeamShapeIds;
00502   double          myPar1[2], myPar2[2]; // U and V bounds of a closed periodic surface
00503   int             myParIndex;     // bounds' index (1-U, 2-V, 3-both)
00504 
00505   typedef std::map< int, GeomAPI_ProjectPointOnSurf* > TID2ProjectorOnSurf;
00506   TID2ProjectorOnSurf myFace2Projector;
00507   typedef std::map< int, GeomAPI_ProjectPointOnCurve* > TID2ProjectorOnCurve;
00508   TID2ProjectorOnCurve myEdge2Projector;
00509 
00510   TopoDS_Shape    myShape;
00511   SMESH_Mesh*     myMesh;
00512   int             myShapeID;
00513 
00514   // to create quadratic elements
00515   bool            myCreateQuadratic;
00516   bool            mySetElemOnShape;
00517 
00518   std::map< int,bool > myNodePosShapesValidity;
00519   bool toCheckPosOnShape(int shapeID ) const;
00520   void setPosOnShapeValidity(int shapeID, bool ok ) const;
00521 
00522 };
00523 
00524 
00525 #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