]> SALOME platform Git repositories - modules/smesh.git/blob - src/SMDS/SMDS_MeshElement.hxx
Salome HOME
Fix
[modules/smesh.git] / src / SMDS / SMDS_MeshElement.hxx
1 // Copyright (C) 2007-2013  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
23 //  SMESH SMDS : implementaion of Salome mesh data structure
24 //  File   : SMDS_MeshElement.hxx
25 //  Module : SMESH
26 //
27 #ifndef _SMDS_MeshElement_HeaderFile
28 #define _SMDS_MeshElement_HeaderFile
29
30 #include "SMESH_SMDS.hxx"
31         
32 #include "SMDSAbs_ElementType.hxx"
33 #include "SMDS_MeshObject.hxx"
34 #include "SMDS_ElemIterator.hxx"
35 #include "SMDS_MeshElementIDFactory.hxx"
36 #include "SMDS_StdIterator.hxx"
37
38 #include <vector>
39 #include <iostream>
40
41 #include <vtkType.h>
42 #include <vtkCellType.h>
43
44 //typedef unsigned short UShortType;
45 typedef short ShortType;
46 typedef int   LongType;
47
48 class SMDS_MeshNode;
49 class SMDS_MeshEdge;
50 class SMDS_MeshFace;
51 class SMDS_Mesh;
52
53 // ============================================================
54 /*!
55  * \brief Base class for elements
56  */
57 // ============================================================
58
59
60 class SMDS_EXPORT SMDS_MeshElement:public SMDS_MeshObject
61 {
62 public:
63
64   SMDS_ElemIteratorPtr nodesIterator() const;
65   SMDS_ElemIteratorPtr edgesIterator() const;
66   SMDS_ElemIteratorPtr facesIterator() const;
67   virtual SMDS_ElemIteratorPtr elementsIterator(SMDSAbs_ElementType type) const;
68   virtual SMDS_ElemIteratorPtr interlacedNodesElemIterator() const;
69
70   virtual SMDS_NodeIteratorPtr nodeIterator() const;
71   virtual SMDS_NodeIteratorPtr interlacedNodesIterator() const;
72   virtual SMDS_NodeIteratorPtr nodesIteratorToUNV() const;
73
74   // std-like iteration on nodes
75   typedef SMDS_StdIterator< const SMDS_MeshNode*, SMDS_ElemIteratorPtr > iterator;
76   iterator begin_nodes() const { return iterator( nodesIterator() ); }
77   iterator end_nodes()   const { return iterator(); }
78
79   virtual int NbNodes() const;
80   virtual int NbEdges() const;
81   virtual int NbFaces() const;
82   inline int GetID() const { return myID; };
83
84   ///Return the type of the current element
85   virtual SMDSAbs_ElementType GetType() const = 0;
86   virtual SMDSAbs_EntityType GetEntityType() const = 0;
87   virtual SMDSAbs_GeometryType GetGeomType() const = 0;
88   virtual vtkIdType GetVtkType() const = 0;
89   virtual bool IsPoly() const { return false; }
90   virtual bool IsQuadratic() const;
91
92   virtual bool IsMediumNode(const SMDS_MeshNode* node) const;
93   virtual int  NbCornerNodes() const;
94
95   friend SMDS_EXPORT std::ostream & operator <<(std::ostream & OS, const SMDS_MeshElement *);
96   friend SMDS_EXPORT bool SMDS_MeshElementIDFactory::BindID(int ID,SMDS_MeshElement* elem);
97   friend class SMDS_Mesh;
98   friend class SMESHDS_Mesh;
99   friend class SMESHDS_SubMesh;
100   friend class SMDS_MeshElementIDFactory;
101
102   // ===========================
103   //  Access to nodes by index
104   // ===========================
105   /*!
106    * \brief Return node by its index
107     * \param ind - node index
108     * \retval const SMDS_MeshNode* - the node
109    */
110   virtual const SMDS_MeshNode* GetNode(const int ind) const;
111
112   /*!
113    * \brief Return node by its index
114     * \param ind - node index
115     * \retval const SMDS_MeshNode* - the node
116    * 
117    * Index is wrapped if it is out of a valid range
118    */
119   const SMDS_MeshNode* GetNodeWrap(const int ind) const { return GetNode( WrappedIndex( ind )); }
120
121   /*!
122    * \brief Return true if index of node is valid (0 <= ind < NbNodes())
123     * \param ind - node index
124     * \retval bool - index check result
125    */
126   virtual bool IsValidIndex(const int ind) const;
127
128   /*!
129    * \brief Return a valid node index, fixing the given one if necessary
130     * \param ind - node index
131     * \retval int - valid node index
132    */
133   int WrappedIndex(const int ind) const {
134     if ( ind < 0 ) return NbNodes() + ind % NbNodes();
135     if ( ind >= NbNodes() ) return ind % NbNodes();
136     return ind;
137   }
138
139   /*!
140    * \brief Check if a node belongs to the element
141     * \param node - the node to check
142     * \retval int - node index within the element, -1 if not found
143    */
144   virtual int GetNodeIndex( const SMDS_MeshNode* node ) const;
145
146   inline ShortType getMeshId()  const { return myMeshId; }
147   inline LongType  getshapeId() const { return myShapeId; }
148   inline int getIdInShape()     const { return myIdInShape; }
149   inline int getVtkId()         const { return myVtkID; }
150
151   /*!
152    * \brief Filters of elements, to be used with SMDS_SetIterator
153    */
154   struct Filter
155   {
156     virtual bool operator()(const SMDS_MeshElement* e) const = 0;
157     ~Filter() {}
158   };
159   struct NonNullFilter: public Filter
160   {
161     bool operator()(const SMDS_MeshElement* e) const { return e; }
162   };
163   struct TypeFilter : public Filter
164   {
165     SMDSAbs_ElementType _type;
166     TypeFilter( SMDSAbs_ElementType t = SMDSAbs_NbElementTypes ):_type(t) {}
167     bool operator()(const SMDS_MeshElement* e) const { return e && e->GetType() == _type; }
168   };
169   struct EntityFilter : public Filter
170   {
171     SMDSAbs_EntityType _type;
172     EntityFilter( SMDSAbs_EntityType t = SMDSEntity_Last ):_type(t) {}
173     bool operator()(const SMDS_MeshElement* e) const { return e && e->GetEntityType() == _type; }
174   };
175   struct GeomFilter : public Filter
176   {
177     SMDSAbs_GeometryType _type;
178     GeomFilter( SMDSAbs_GeometryType t = SMDSGeom_NONE ):_type(t) {}
179     bool operator()(const SMDS_MeshElement* e) const { return e && e->GetGeomType() == _type; }
180   };
181
182 protected:
183   inline void setId(int id) {myID = id; }
184   inline void setShapeId(LongType shapeId) {myShapeId = shapeId; }
185   inline void setIdInShape(int id) { myIdInShape = id; }
186   inline void setVtkId(int vtkId) { myVtkID = vtkId; }
187   SMDS_MeshElement(int ID=-1);
188   SMDS_MeshElement(int id, ShortType meshId, LongType shapeId = 0);
189   virtual void init(int id = -1, ShortType meshId = -1, LongType shapeId = 0);
190   virtual void Print(std::ostream & OS) const;
191
192   //! Element index in vector SMDS_Mesh::myNodes or SMDS_Mesh::myCells
193   int myID;
194   //! index in vtkUnstructuredGrid
195   int myVtkID;
196   //! SMDS_Mesh identification in SMESH
197   ShortType myMeshId;
198   //! SubShape and SubMesh identification in SMESHDS
199   LongType myShapeId;
200   //! Element index in SMESHDS_SubMesh vector
201   int myIdInShape;
202 };
203
204
205 // ============================================================
206 /*!
207  * \brief Comparator of elements by ID for usage in std containers
208  */
209 // ============================================================
210
211 struct TIDTypeCompare {
212   bool operator () (const SMDS_MeshElement* e1, const SMDS_MeshElement* e2) const
213   { return e1->GetType() == e2->GetType() ? e1->GetID() < e2->GetID() : e1->GetType() < e2->GetType(); }
214 };
215
216 // WARNING: this comparator makes impossible to store both nodes and elements in the same set
217 // because there are nodes and elements with the same ID. Use TIDTypeCompare for such containers.
218 struct TIDCompare {
219   bool operator () (const SMDS_MeshElement* e1, const SMDS_MeshElement* e2) const
220   { return e1->GetID() < e2->GetID(); }
221 };
222
223 #endif