Back to index

salome-smesh  6.5.0
SMESH_MesherHelper.hxx
Go to the documentation of this file.
00001 // Copyright (C) 2007-2012  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 std::list<TopoDS_Edge>& theBaseSide,
00100                               SMESHDS_Mesh*                 theMesh,
00101                               SMESH_ProxyMesh*              theProxyMesh=0);
00105   static bool LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
00106                               const TopoDS_Face& theFace,
00107                               const TopoDS_Edge& theBaseEdge,
00108                               SMESHDS_Mesh*      theMesh,
00109                               SMESH_ProxyMesh*   theProxyMesh=0);
00116   static TopoDS_Shape GetSubShapeByNode(const SMDS_MeshNode* node,
00117                                         const SMESHDS_Mesh*  meshDS);
00118 
00125   static int WrapIndex(const int ind, const int nbNodes) {
00126     if ( ind < 0 ) return nbNodes + ind % nbNodes;
00127     if ( ind >= nbNodes ) return ind % nbNodes;
00128     return ind;
00129   }
00130 
00134   static int NbAncestors(const TopoDS_Shape& shape,
00135                          const SMESH_Mesh&   mesh,
00136                          TopAbs_ShapeEnum    ancestorType=TopAbs_SHAPE);
00140   static PShapeIteratorPtr GetAncestors(const TopoDS_Shape& shape,
00141                                         const SMESH_Mesh&   mesh,
00142                                         TopAbs_ShapeEnum    ancestorType);
00146   static TopoDS_Shape GetCommonAncestor(const TopoDS_Shape& shape1,
00147                                         const TopoDS_Shape& shape2,
00148                                         const SMESH_Mesh&   mesh,
00149                                         TopAbs_ShapeEnum    ancestorType);
00150 
00154   static TopAbs_Orientation GetSubShapeOri(const TopoDS_Shape& shape,
00155                                            const TopoDS_Shape& subShape);
00156 
00157   static bool IsSubShape( const TopoDS_Shape& shape, const TopoDS_Shape& mainShape );
00158 
00159   static bool IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh );
00160 
00161   static double MaxTolerance( const TopoDS_Shape& shape );
00162 
00163   static bool IsClosedEdge( const TopoDS_Edge& anEdge );
00164 
00165   static TopoDS_Vertex IthVertex( const bool is2nd, TopoDS_Edge anEdge, const bool CumOri=true );
00166 
00167 
00168 public:
00169   // ---------- PUBLIC INSTANCE METHODS ----------
00170 
00171   // constructor
00172   SMESH_MesherHelper(SMESH_Mesh& theMesh);
00173 
00174   SMESH_Mesh* GetMesh() const { return myMesh; }
00175     
00176   SMESHDS_Mesh* GetMeshDS() const { return GetMesh()->GetMeshDS(); }
00177     
00182   bool IsQuadraticSubMesh(const TopoDS_Shape& theShape);
00186   void SetIsQuadratic(const bool theBuildQuadratic)
00187   { myCreateQuadratic = theBuildQuadratic; }
00191   bool GetIsQuadratic() const { return myCreateQuadratic; }
00192 
00197   void FixQuadraticElements(bool volumeOnly=true);
00198 
00204   void SetElementsOnShape(bool toSet) { mySetElemOnShape = toSet; }
00205 
00209   void SetSubShape(const int           subShapeID);
00210   void SetSubShape(const TopoDS_Shape& subShape);
00215   int GetSubShapeID() const { return myShapeID; }
00219   const TopoDS_Shape& GetSubShape() const  { return myShape; }
00220 
00224   SMDS_MeshNode* AddNode(double x, double y, double z, int ID = 0);
00228   SMDS_MeshEdge* AddEdge(const SMDS_MeshNode* n1,
00229                          const SMDS_MeshNode* n2,
00230                          const int id = 0, 
00231                          const bool force3d = true);
00235   SMDS_MeshFace* AddFace(const SMDS_MeshNode* n1,
00236                          const SMDS_MeshNode* n2,
00237                          const SMDS_MeshNode* n3,
00238                          const int id=0, 
00239                          const bool force3d = false);
00243   SMDS_MeshFace* AddFace(const SMDS_MeshNode* n1,
00244                          const SMDS_MeshNode* n2,
00245                          const SMDS_MeshNode* n3,
00246                          const SMDS_MeshNode* n4,
00247                          const int id = 0,
00248                          const bool force3d = false);
00249 
00253   SMDS_MeshFace* AddPolygonalFace (const std::vector<const SMDS_MeshNode*>& nodes,
00254                                    const int id = 0,
00255                                    const bool force3d = false);
00259   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00260                              const SMDS_MeshNode* n2,
00261                              const SMDS_MeshNode* n3,
00262                              const SMDS_MeshNode* n4,
00263                              const int id = 0,
00264                              const bool force3d = true);
00268   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00269                              const SMDS_MeshNode* n2,
00270                              const SMDS_MeshNode* n3,
00271                              const SMDS_MeshNode* n4,
00272                              const SMDS_MeshNode* n5,
00273                              const int id = 0,
00274                              const bool force3d = true);
00278   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00279                              const SMDS_MeshNode* n2,
00280                              const SMDS_MeshNode* n3,
00281                              const SMDS_MeshNode* n4,
00282                              const SMDS_MeshNode* n5,
00283                              const SMDS_MeshNode* n6,
00284                              const int id = 0, 
00285                              const bool force3d = true);
00289   SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
00290                              const SMDS_MeshNode* n2,
00291                              const SMDS_MeshNode* n3,
00292                              const SMDS_MeshNode* n4,
00293                              const SMDS_MeshNode* n5,
00294                              const SMDS_MeshNode* n6,
00295                              const SMDS_MeshNode* n7,
00296                              const SMDS_MeshNode* n8,
00297                              const int id = 0, 
00298                              bool force3d = true);
00299 
00303   SMDS_MeshVolume* AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
00304                                         const std::vector<int>&                  quantities,
00305                                         const int                                ID=0,
00306                                         const bool                               force3d = true);
00310   double GetNodeU(const TopoDS_Edge&   theEdge,
00311                   const SMDS_MeshNode* theNode,
00312                   const SMDS_MeshNode* inEdgeNode=0,
00313                   bool*                check=0);
00318   gp_XY GetNodeUV(const TopoDS_Face&   F,
00319                   const SMDS_MeshNode* n,
00320                   const SMDS_MeshNode* inFaceNode=0,
00321                   bool*                check=0) const;
00328   bool CheckNodeUV(const TopoDS_Face&   F,
00329                    const SMDS_MeshNode* n,
00330                    gp_XY&               uv,
00331                    const double         tol,
00332                    const bool           force=false,
00333                    double               distXYZ[4]=0) const;
00340   bool CheckNodeU(const TopoDS_Edge&   E,
00341                   const SMDS_MeshNode* n,
00342                   double&              u,
00343                   const double         tol,
00344                   const bool           force=false,
00345                   double               distXYZ[4]=0) const;
00349   static gp_XY GetMiddleUV(const Handle(Geom_Surface)& surface,
00350                            const gp_XY&                uv1,
00351                            const gp_XY&                uv2);
00359 #define gp_XY_FunPtr(meth) \
00360   static gp_XY __gpXY_##meth (const gp_XY& uv1, const gp_XY& uv2) { return uv1.meth( uv2 ); } \
00361   static xyFunPtr gp_XY_##meth = & __gpXY_##meth
00362 
00368   static gp_XY applyIn2D(const Handle(Geom_Surface)& surface,
00369                          const gp_XY&                uv1,
00370                          const gp_XY&                uv2,
00371                          xyFunPtr                    fun,
00372                          const bool                  resultInPeriod=true);
00373                           
00381   bool GetNodeUVneedInFaceNode(const TopoDS_Face& F = TopoDS_Face()) const;
00382 
00386   GeomAPI_ProjectPointOnSurf& GetProjector(const TopoDS_Face& F,
00387                                            TopLoc_Location&   loc,
00388                                            double             tol=0 ) const; 
00389 
00397   bool IsDegenShape(const int subShape) const
00398   { return myDegenShapeIds.find( subShape ) != myDegenShapeIds.end(); }
00404   bool HasDegeneratedEdges() const { return !myDegenShapeIds.empty(); }
00405 
00414   bool IsSeamShape(const int subShape) const
00415   { return mySeamShapeIds.find( subShape ) != mySeamShapeIds.end(); }
00424   bool IsSeamShape(const TopoDS_Shape& subShape) const
00425   { return IsSeamShape( GetMeshDS()->ShapeToIndex( subShape )); }
00430   bool IsRealSeam(const int subShape) const
00431   { return mySeamShapeIds.find( -subShape ) != mySeamShapeIds.end(); }
00436   bool IsRealSeam(const TopoDS_Shape& subShape) const
00437   { return IsRealSeam( GetMeshDS()->ShapeToIndex( subShape)); }
00443   bool HasSeam() const { return !mySeamShapeIds.empty(); }
00448   int GetPeriodicIndex() const { return myParIndex; }
00452   double GetOtherParam(const double param) const;
00453 
00460   const SMDS_MeshNode* GetMediumNode(const SMDS_MeshNode* n1,
00461                                      const SMDS_MeshNode* n2,
00462                                      const bool force3d);
00466   std::pair<int, TopAbs_ShapeEnum> GetMediumPos(const SMDS_MeshNode* n1,
00467                                                 const SMDS_MeshNode* n2);
00471   void AddTLinkNode(const SMDS_MeshNode* n1,
00472                     const SMDS_MeshNode* n2,
00473                     const SMDS_MeshNode* n12);
00477   void AddTLinkNodeMap(const TLinkNodeMap& aMap)
00478     { myTLinkNodeMap.insert(aMap.begin(), aMap.end()); }
00479 
00480   void AddTLinks(const SMDS_MeshEdge*   edge);
00481   void AddTLinks(const SMDS_MeshFace*   face);
00482   void AddTLinks(const SMDS_MeshVolume* vol);
00483 
00487   const TLinkNodeMap& GetTLinkNodeMap() const { return myTLinkNodeMap; }
00488 
00494   enum MType{ LINEAR, QUADRATIC, COMP };
00495   MType IsQuadraticMesh();
00496   
00497   virtual ~SMESH_MesherHelper();
00498 
00499 protected:
00500 
00507   gp_Pnt2d GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const;
00508 
00509   const SMDS_MeshNode* getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
00510                                                    const SMDS_MeshNode* n2,
00511                                                    bool                 force3d);
00512  private:
00513 
00514   // Forbiden copy constructor
00515   SMESH_MesherHelper (const SMESH_MesherHelper& theOther) {};
00516 
00517   // special map for using during creation of quadratic elements
00518   TLinkNodeMap    myTLinkNodeMap;
00519 
00520   std::set< int > myDegenShapeIds;
00521   std::set< int > mySeamShapeIds;
00522   double          myPar1[2], myPar2[2]; // U and V bounds of a closed periodic surface
00523   int             myParIndex;     // bounds' index (1-U, 2-V, 3-both)
00524 
00525   typedef std::map< int, GeomAPI_ProjectPointOnSurf* > TID2ProjectorOnSurf;
00526   TID2ProjectorOnSurf myFace2Projector;
00527   typedef std::map< int, GeomAPI_ProjectPointOnCurve* > TID2ProjectorOnCurve;
00528   TID2ProjectorOnCurve myEdge2Projector;
00529 
00530   TopoDS_Shape    myShape;
00531   SMESH_Mesh*     myMesh;
00532   int             myShapeID;
00533 
00534   // to create quadratic elements
00535   bool            myCreateQuadratic;
00536   bool            mySetElemOnShape;
00537 
00538   std::map< int,bool > myNodePosShapesValidity;
00539   bool toCheckPosOnShape(int shapeID ) const;
00540   void setPosOnShapeValidity(int shapeID, bool ok ) const;
00541 
00542 };
00543 
00544 
00545 #endif