Salome HOME
PR: merged from V5_1_4rc1
[modules/smesh.git] / src / SMDS / SMDS_MeshElementIDFactory.cxx
1 //  Copyright (C) 2007-2010  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_MeshElementIDFactory.cxx
25 //  Author : Jean-Michel BOULCOURT
26 //  Module : SMESH
27 //
28 #ifdef _MSC_VER
29 #pragma warning(disable:4786)
30 #endif
31
32 #include "SMDS_MeshElementIDFactory.hxx"
33 #include "SMDS_MeshElement.hxx"
34 #include "SMDS_Mesh.hxx"
35
36 #include "utilities.h"
37
38 #include <vtkUnstructuredGrid.h>
39 #include <vtkCellType.h>
40
41 using namespace std;
42
43 //=======================================================================
44 //function : SMDS_MeshElementIDFactory
45 //purpose  : 
46 //=======================================================================
47 SMDS_MeshElementIDFactory::SMDS_MeshElementIDFactory():
48   SMDS_MeshNodeIDFactory()
49 {
50 //    myIDElements.clear();
51 //    myVtkIndex.clear();
52     myVtkCellTypes.clear();
53     myVtkCellTypes.reserve(SMDSEntity_Last);
54     myVtkCellTypes[SMDSEntity_Node]            = VTK_VERTEX;
55     myVtkCellTypes[SMDSEntity_0D]              = VTK_VERTEX;
56     myVtkCellTypes[SMDSEntity_Edge]            = VTK_LINE;
57     myVtkCellTypes[SMDSEntity_Quad_Edge]       = VTK_QUADRATIC_EDGE;
58     myVtkCellTypes[SMDSEntity_Triangle]        = VTK_TRIANGLE;
59     myVtkCellTypes[SMDSEntity_Quad_Triangle]   = VTK_QUADRATIC_TRIANGLE;
60     myVtkCellTypes[SMDSEntity_Quadrangle]      = VTK_QUAD;
61     myVtkCellTypes[SMDSEntity_Quad_Quadrangle] = VTK_QUADRATIC_TRIANGLE;
62     myVtkCellTypes[SMDSEntity_Polygon]         = VTK_POLYGON;
63     myVtkCellTypes[SMDSEntity_Quad_Polygon]    = VTK_POLYGON; // -PR- verifer
64     myVtkCellTypes[SMDSEntity_Tetra]           = VTK_TETRA;
65     myVtkCellTypes[SMDSEntity_Quad_Tetra]      = VTK_QUADRATIC_TETRA;
66     myVtkCellTypes[SMDSEntity_Pyramid]         = VTK_PYRAMID;
67     myVtkCellTypes[SMDSEntity_Quad_Pyramid]    = VTK_CONVEX_POINT_SET;
68     myVtkCellTypes[SMDSEntity_Hexa]            = VTK_HEXAHEDRON;
69     myVtkCellTypes[SMDSEntity_Quad_Hexa]       = VTK_QUADRATIC_HEXAHEDRON;
70     myVtkCellTypes[SMDSEntity_Penta]           = VTK_WEDGE;
71     myVtkCellTypes[SMDSEntity_Quad_Penta]      = VTK_QUADRATIC_WEDGE;
72     myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_CONVEX_POINT_SET;
73     myVtkCellTypes[SMDSEntity_Quad_Polyhedra]  = VTK_CONVEX_POINT_SET;
74 }
75
76 int SMDS_MeshElementIDFactory::SetInVtkGrid(SMDS_MeshElement * elem)
77 {
78    // --- retreive nodes ID
79
80   vector<vtkIdType> nodeIds;
81   SMDS_ElemIteratorPtr it = elem->nodesIterator();
82   while(it->more())
83   {
84       int nodeId = it->next()->GetID();
85       //MESSAGE("   node in cell " << ID << " : " <<nodeId)
86       nodeIds.push_back(nodeId);
87   }
88
89   // --- insert cell in vtkUnstructuredGrid
90
91   vtkUnstructuredGrid * grid = myMesh->getGrid();
92   int typ = GetVtkCellType(elem->GetType());
93   int cellId = grid->InsertNextLinkedCell(typ, nodeIds.size(), &nodeIds[0]);
94   SMDS_MeshCell *cell = dynamic_cast<SMDS_MeshCell*>(elem);
95   assert(cell);
96   cell->setVtkId(cellId); 
97   //MESSAGE("SMDS_MeshElementIDFactory::SetInVtkGrid " << cellId);
98   return cellId;
99 }
100
101 //=======================================================================
102 //function : BindID
103 //purpose  : 
104 //=======================================================================
105
106 bool SMDS_MeshElementIDFactory::BindID(int ID, SMDS_MeshElement * elem)
107 {
108   MESSAGE("SMDS_MeshElementIDFactory::BindID " << ID);
109   SetInVtkGrid(elem);
110   return myMesh->registerElement(ID, elem);
111 }
112
113 //=======================================================================
114 //function : MeshElement
115 //purpose  : 
116 //=======================================================================
117 SMDS_MeshElement* SMDS_MeshElementIDFactory::MeshElement(int ID)
118 {
119   if ((ID<0) || (ID>myMax) || (myMesh->myIDElements[ID]<0))
120     return NULL;
121   const SMDS_MeshElement* elem = GetMesh()->FindElement(ID);
122   return (SMDS_MeshElement*)(elem);
123 }
124
125 //=======================================================================
126 //function : ReleaseID
127 //purpose  : 
128 //=======================================================================
129 void SMDS_MeshElementIDFactory::ReleaseID(const int ID)
130 {
131   int vtkId = myMesh->myIDElements[ID];
132   //MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID smdsId vtkId " << ID << " " << vtkId);
133   myMesh->myIDElements[ID] = -1;
134   myMesh->myVtkIndex[vtkId] = -1;
135   SMDS_MeshIDFactory::ReleaseID(ID);
136   if (ID == myMax)
137     myMax = 0;
138   if (ID == myMin)
139     myMax = 0;
140 }
141
142 //=======================================================================
143 //function : updateMinMax
144 //purpose  : 
145 //=======================================================================
146
147 void SMDS_MeshElementIDFactory::updateMinMax() const
148 {
149   myMin = IntegerLast();
150   myMax = 0;
151   for (int i=0; i<myMesh->myIDElements.size(); i++)
152       if (int id=myMesh->myIDElements[i] >=0)
153       {
154         if (id > myMax) myMax = id;
155         if (id < myMin) myMin = id;
156       }
157   if (myMin == IntegerLast())
158     myMin = 0;
159 }
160
161 //=======================================================================
162 //function : elementsIterator
163 //purpose  : Return an iterator on elements of the factory
164 //=======================================================================
165
166 SMDS_ElemIteratorPtr SMDS_MeshElementIDFactory::elementsIterator() const
167 {
168     return myMesh->elementsIterator(SMDSAbs_All);
169 }
170
171 void SMDS_MeshElementIDFactory::Clear()
172 {
173   myMesh->myIDElements.clear();
174   myMin = myMax = 0;
175   SMDS_MeshIDFactory::Clear();
176 }
177
178 int SMDS_MeshElementIDFactory::GetVtkCellType(int SMDSType)
179 {
180     assert((SMDSType >=0) && (SMDSType< SMDSEntity_Last));
181     return myVtkCellTypes[SMDSType];
182 }