Salome HOME
1abf9e81f24928a9d0c42c65c1d3cdc6f624d1c5
[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 "SMDS_UnstructuredGrid.hxx"
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 //#ifdef VTK_HAVE_POLYHEDRON
73     myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_POLYHEDRON;
74 //#else
75 //    myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_CONVEX_POINT_SET;
76 //#endif
77     myVtkCellTypes[SMDSEntity_Quad_Polyhedra]  = VTK_CONVEX_POINT_SET;
78 }
79
80 int SMDS_MeshElementIDFactory::SetInVtkGrid(SMDS_MeshElement * elem)
81 {
82    // --- retrieve nodes ID
83
84   SMDS_MeshCell *cell = dynamic_cast<SMDS_MeshCell*>(elem);
85   assert(cell);
86   vector<vtkIdType> nodeIds;
87   SMDS_ElemIteratorPtr it = elem->nodesIterator();
88   while(it->more())
89   {
90       int nodeId = (static_cast<const SMDS_MeshNode*>(it->next()))->getVtkId();
91       MESSAGE("   node in cell " << cell->getVtkId() << " : " << nodeId)
92       nodeIds.push_back(nodeId);
93   }
94
95   // --- insert cell in vtkUnstructuredGrid
96
97   vtkUnstructuredGrid * grid = myMesh->getGrid();
98   int typ = GetVtkCellType(elem->GetType());
99   int cellId = grid->InsertNextLinkedCell(typ, nodeIds.size(), &nodeIds[0]);
100   cell->setVtkId(cellId); 
101   //MESSAGE("SMDS_MeshElementIDFactory::SetInVtkGrid " << cellId);
102   return cellId;
103 }
104
105 //=======================================================================
106 //function : BindID
107 //purpose  : 
108 //=======================================================================
109
110 bool SMDS_MeshElementIDFactory::BindID(int ID, SMDS_MeshElement * elem)
111 {
112   MESSAGE("SMDS_MeshElementIDFactory::BindID " << ID);
113   SetInVtkGrid(elem);
114   return myMesh->registerElement(ID, elem);
115 }
116
117 //=======================================================================
118 //function : MeshElement
119 //purpose  : 
120 //=======================================================================
121 SMDS_MeshElement* SMDS_MeshElementIDFactory::MeshElement(int ID)
122 {
123   if ((ID<1) || (ID>=myMesh->myCells.size()))
124     return NULL;
125   const SMDS_MeshElement* elem = GetMesh()->FindElement(ID);
126   return (SMDS_MeshElement*)(elem);
127 }
128
129 //=======================================================================
130 //function : GetFreeID
131 //purpose  : 
132 //=======================================================================
133
134 int SMDS_MeshElementIDFactory::GetFreeID()
135 {
136   int ID;
137   do {
138     ID = SMDS_MeshIDFactory::GetFreeID();
139   } while ( MeshElement( ID ));
140   return ID;
141 }
142
143 //=======================================================================
144 //function : ReleaseID
145 //purpose  : 
146 //=======================================================================
147 void SMDS_MeshElementIDFactory::ReleaseID(int ID, int vtkId)
148 {
149   if (ID < 1) // TODO check case ID == O
150     {
151       MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID ID = " << ID);
152       return;
153     }
154   //MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID smdsId vtkId " << ID << " " << vtkId);
155   if (vtkId >= 0)
156     {
157       assert(vtkId < myMesh->myCellIdVtkToSmds.size());
158       myMesh->myCellIdVtkToSmds[vtkId] = -1;
159       myMesh->setMyModified();
160     }
161   SMDS_MeshIDFactory::ReleaseID(ID);
162   if (ID == myMax)
163     myMax = 0;
164   if (ID == myMin)
165     myMax = 0;
166 }
167
168 //=======================================================================
169 //function : updateMinMax
170 //purpose  : 
171 //=======================================================================
172
173 void SMDS_MeshElementIDFactory::updateMinMax() const
174 {
175   myMin = IntegerLast();
176   myMax = 0;
177   for (int i = 0; i < myMesh->myCells.size(); i++)
178     {
179       if (myMesh->myCells[i])
180         {
181           int id = myMesh->myCells[i]->GetID();
182           if (id > myMax)
183             myMax = id;
184           if (id < myMin)
185             myMin = id;
186         }
187     }
188   if (myMin == IntegerLast())
189     myMin = 0;
190 }
191
192 //=======================================================================
193 //function : elementsIterator
194 //purpose  : Return an iterator on elements of the factory
195 //=======================================================================
196
197 SMDS_ElemIteratorPtr SMDS_MeshElementIDFactory::elementsIterator() const
198 {
199     return myMesh->elementsIterator(SMDSAbs_All);
200 }
201
202 void SMDS_MeshElementIDFactory::Clear()
203 {
204   //myMesh->myCellIdSmdsToVtk.clear();
205   myMesh->myCellIdVtkToSmds.clear();
206   myMin = myMax = 0;
207   SMDS_MeshIDFactory::Clear();
208 }
209
210 int SMDS_MeshElementIDFactory::GetVtkCellType(int SMDSType)
211 {
212     assert((SMDSType >=0) && (SMDSType< SMDSEntity_Last));
213     return myVtkCellTypes[SMDSType];
214 }