Salome HOME
Update French translation file
[modules/smesh.git] / src / SMDS / SMDS_MeshElementIDFactory.cxx
1 // Copyright (C) 2007-2011  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 #include <climits>
42
43 using namespace std;
44
45 //=======================================================================
46 //function : SMDS_MeshElementIDFactory
47 //purpose  : 
48 //=======================================================================
49 SMDS_MeshElementIDFactory::SMDS_MeshElementIDFactory():
50   SMDS_MeshNodeIDFactory()
51 {
52 //    myIDElements.clear();
53 //    myVtkIndex.clear();
54     myVtkCellTypes.clear();
55     myVtkCellTypes.reserve(SMDSEntity_Last);
56     myVtkCellTypes[SMDSEntity_Node]            = VTK_VERTEX;
57     myVtkCellTypes[SMDSEntity_0D]              = VTK_VERTEX;
58     myVtkCellTypes[SMDSEntity_Edge]            = VTK_LINE;
59     myVtkCellTypes[SMDSEntity_Quad_Edge]       = VTK_QUADRATIC_EDGE;
60     myVtkCellTypes[SMDSEntity_Triangle]        = VTK_TRIANGLE;
61     myVtkCellTypes[SMDSEntity_Quad_Triangle]   = VTK_QUADRATIC_TRIANGLE;
62     myVtkCellTypes[SMDSEntity_Quadrangle]      = VTK_QUAD;
63     myVtkCellTypes[SMDSEntity_Quad_Quadrangle] = VTK_QUADRATIC_TRIANGLE;
64     myVtkCellTypes[SMDSEntity_Polygon]         = VTK_POLYGON;
65     myVtkCellTypes[SMDSEntity_Quad_Polygon]    = VTK_POLYGON; // -PR- verifer
66     myVtkCellTypes[SMDSEntity_Tetra]           = VTK_TETRA;
67     myVtkCellTypes[SMDSEntity_Quad_Tetra]      = VTK_QUADRATIC_TETRA;
68     myVtkCellTypes[SMDSEntity_Pyramid]         = VTK_PYRAMID;
69     myVtkCellTypes[SMDSEntity_Quad_Pyramid]    = VTK_CONVEX_POINT_SET;
70     myVtkCellTypes[SMDSEntity_Hexa]            = VTK_HEXAHEDRON;
71     myVtkCellTypes[SMDSEntity_Quad_Hexa]       = VTK_QUADRATIC_HEXAHEDRON;
72     myVtkCellTypes[SMDSEntity_Penta]           = VTK_WEDGE;
73     myVtkCellTypes[SMDSEntity_Quad_Penta]      = VTK_QUADRATIC_WEDGE;
74 //#ifdef VTK_HAVE_POLYHEDRON
75     myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_POLYHEDRON;
76 //#else
77 //    myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_CONVEX_POINT_SET;
78 //#endif
79     myVtkCellTypes[SMDSEntity_Quad_Polyhedra]  = VTK_CONVEX_POINT_SET;
80 }
81
82 int SMDS_MeshElementIDFactory::SetInVtkGrid(SMDS_MeshElement * elem)
83 {
84    // --- retrieve nodes ID
85
86   SMDS_MeshCell *cell = dynamic_cast<SMDS_MeshCell*>(elem);
87   assert(cell);
88   vector<vtkIdType> nodeIds;
89   SMDS_ElemIteratorPtr it = elem->nodesIterator();
90   while(it->more())
91   {
92       int nodeId = (static_cast<const SMDS_MeshNode*>(it->next()))->getVtkId();
93       MESSAGE("   node in cell " << cell->getVtkId() << " : " << nodeId)
94       nodeIds.push_back(nodeId);
95   }
96
97   // --- insert cell in vtkUnstructuredGrid
98
99   vtkUnstructuredGrid * grid = myMesh->getGrid();
100   //int locType = elem->GetType();
101   int typ = VTK_VERTEX;//GetVtkCellType(locType);
102   int cellId = grid->InsertNextLinkedCell(typ, nodeIds.size(), &nodeIds[0]);
103   cell->setVtkId(cellId); 
104   //MESSAGE("SMDS_MeshElementIDFactory::SetInVtkGrid " << cellId);
105   return cellId;
106 }
107
108 //=======================================================================
109 //function : BindID
110 //purpose  : 
111 //=======================================================================
112
113 bool SMDS_MeshElementIDFactory::BindID(int ID, SMDS_MeshElement * elem)
114 {
115   MESSAGE("SMDS_MeshElementIDFactory::BindID " << ID);
116   SetInVtkGrid(elem);
117   return myMesh->registerElement(ID, elem);
118 }
119
120 //=======================================================================
121 //function : MeshElement
122 //purpose  : 
123 //=======================================================================
124 SMDS_MeshElement* SMDS_MeshElementIDFactory::MeshElement(int ID)
125 {
126   if ((ID<1) || (ID>=myMesh->myCells.size()))
127     return NULL;
128   const SMDS_MeshElement* elem = GetMesh()->FindElement(ID);
129   return (SMDS_MeshElement*)(elem);
130 }
131
132 //=======================================================================
133 //function : GetFreeID
134 //purpose  : 
135 //=======================================================================
136
137 int SMDS_MeshElementIDFactory::GetFreeID()
138 {
139   int ID;
140   do {
141     ID = SMDS_MeshIDFactory::GetFreeID();
142   } while ( MeshElement( ID ));
143   return ID;
144 }
145
146 //=======================================================================
147 //function : ReleaseID
148 //purpose  : 
149 //=======================================================================
150 void SMDS_MeshElementIDFactory::ReleaseID(int ID, int vtkId)
151 {
152   if (ID < 1) // TODO check case ID == O
153     {
154       MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID ID = " << ID);
155       return;
156     }
157   //MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID smdsId vtkId " << ID << " " << vtkId);
158   if (vtkId >= 0)
159     {
160       assert(vtkId < myMesh->myCellIdVtkToSmds.size());
161       myMesh->myCellIdVtkToSmds[vtkId] = -1;
162       myMesh->setMyModified();
163     }
164   SMDS_MeshIDFactory::ReleaseID(ID);
165   if (ID == myMax)
166     myMax = 0;
167   if (ID == myMin)
168     myMax = 0;
169 }
170
171 //=======================================================================
172 //function : updateMinMax
173 //purpose  : 
174 //=======================================================================
175
176 void SMDS_MeshElementIDFactory::updateMinMax() const
177 {
178   myMin = INT_MAX;
179   myMax = 0;
180   for (int i = 0; i < myMesh->myCells.size(); i++)
181     {
182       if (myMesh->myCells[i])
183         {
184           int id = myMesh->myCells[i]->GetID();
185           if (id > myMax)
186             myMax = id;
187           if (id < myMin)
188             myMin = id;
189         }
190     }
191   if (myMin == INT_MAX)
192     myMin = 0;
193 }
194
195 //=======================================================================
196 //function : elementsIterator
197 //purpose  : Return an iterator on elements of the factory
198 //=======================================================================
199
200 SMDS_ElemIteratorPtr SMDS_MeshElementIDFactory::elementsIterator() const
201 {
202     return myMesh->elementsIterator(SMDSAbs_All);
203 }
204
205 void SMDS_MeshElementIDFactory::Clear()
206 {
207   //myMesh->myCellIdSmdsToVtk.clear();
208   myMesh->myCellIdVtkToSmds.clear();
209   myMin = myMax = 0;
210   SMDS_MeshIDFactory::Clear();
211 }
212
213 int SMDS_MeshElementIDFactory::GetVtkCellType(int SMDSType)
214 {
215     assert((SMDSType >=0) && (SMDSType< SMDSEntity_Last));
216     return myVtkCellTypes[SMDSType];
217 }