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