Salome HOME
23315: [CEA 1929] Too much memory used to display a mesh in shading and wireframe
[modules/smesh.git] / src / SMDS / SMDS_MeshNode.cxx
1 // Copyright (C) 2007-2016  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, or (at your option) any later version.
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 //
25 #ifdef _MSC_VER
26 #pragma warning(disable:4786)
27 #endif
28
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_SpacePosition.hxx"
31 #include "SMDS_IteratorOfElements.hxx"
32 #include "SMDS_Mesh.hxx"
33 #include <vtkUnstructuredGrid.h>
34
35 #include "utilities.h"
36 #include "Utils_SALOME_Exception.hxx"
37 #include <cassert>
38
39 using namespace std;
40
41 int SMDS_MeshNode::nbNodes =0;
42
43 //=======================================================================
44 //function : SMDS_MeshNode
45 //purpose  :
46 //=======================================================================
47 SMDS_MeshNode::SMDS_MeshNode() :
48   SMDS_MeshElement(-1, -1, 0),
49   myPosition(SMDS_SpacePosition::originSpacePosition())
50 {
51   nbNodes++;
52 }
53
54 SMDS_MeshNode::SMDS_MeshNode(int id, int meshId, int shapeId, double x, double y, double z):
55   SMDS_MeshElement(id, meshId, shapeId),
56   myPosition(SMDS_SpacePosition::originSpacePosition())
57 {
58   nbNodes++;
59   init(id, meshId, shapeId, x, y ,z);
60 }
61
62 void SMDS_MeshNode::init(int id, int meshId, int shapeId, double x, double y, double z)
63 {
64   SMDS_MeshElement::init(id, meshId, shapeId);
65   myVtkID = id - 1;
66   assert(myVtkID >= 0);
67   SMDS_UnstructuredGrid * grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
68   vtkPoints *points = grid->GetPoints();
69   points->InsertPoint(myVtkID, x, y, z);
70   grid->GetLinks()->ResizeForPoint( myVtkID );
71 }
72
73 SMDS_MeshNode::~SMDS_MeshNode()
74 {
75   nbNodes--;
76   if ( myPosition && myPosition != SMDS_SpacePosition::originSpacePosition() )
77     delete myPosition, myPosition = 0;
78 }
79
80 //=======================================================================
81 //function : RemoveInverseElement
82 //purpose  :
83 //=======================================================================
84
85 void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * parent)
86 {
87   //MESSAGE("RemoveInverseElement " << myID << " " << parent->GetID());
88   const SMDS_MeshCell* cell = dynamic_cast<const SMDS_MeshCell*>(parent);
89   MYASSERT(cell);
90   if ( SMDS_Mesh::_meshList[myMeshId]->getGrid()->HasLinks() )
91     SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, cell->getVtkId());
92 }
93
94 //=======================================================================
95 //function : Print
96 //purpose  :
97 //=======================================================================
98
99 void SMDS_MeshNode::Print(ostream & OS) const
100 {
101   OS << "Node <" << myID << "> : X = " << X() << " Y = "
102      << Y() << " Z = " << Z() << endl;
103 }
104
105 //=======================================================================
106 //function : SetPosition
107 //purpose  :
108 //=======================================================================
109
110 void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
111 {
112   if ( myPosition &&
113        myPosition != SMDS_SpacePosition::originSpacePosition() &&
114        myPosition != aPos )
115     delete myPosition;
116   myPosition = aPos;
117 }
118
119 //=======================================================================
120 //function : GetPosition
121 //purpose  :
122 //=======================================================================
123
124 const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
125 {
126   return myPosition;
127 }
128
129 //=======================================================================
130 /*!
131  * \brief Iterator on list of elements
132  */
133 //=======================================================================
134
135 class SMDS_MeshNode_MyInvIterator: public SMDS_ElemIterator
136 {
137 private:
138   SMDS_Mesh* myMesh;
139   vtkIdType* myCells;
140   int myNcells;
141   SMDSAbs_ElementType myType;
142   int iter;
143   vector<vtkIdType> cellList;
144
145 public:
146   SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
147     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
148   {
149     cellList.reserve( ncells );
150     if (type == SMDSAbs_All)
151       cellList.assign( cells, cells + ncells );
152     else
153       for (int i = 0; i < ncells; i++)
154       {
155         int vtkId = cells[i];
156         int smdsId = myMesh->fromVtkToSmds(vtkId);
157         const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
158         if (elem->GetType() == type)
159         {
160           cellList.push_back(vtkId);
161         }
162       }
163     myCells = cellList.empty() ? 0 : &cellList[0];
164     myNcells = cellList.size();
165   }
166
167   bool more()
168   {
169     return (iter < myNcells);
170   }
171
172   const SMDS_MeshElement* next()
173   {
174     int vtkId = myCells[iter];
175     int smdsId = myMesh->fromVtkToSmds(vtkId);
176     const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
177     if (!elem)
178     {
179       MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
180       throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
181     }
182     iter++;
183     return elem;
184   }
185 };
186
187 SMDS_ElemIteratorPtr SMDS_MeshNode::GetInverseElementIterator(SMDSAbs_ElementType type) const
188 {
189   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
190   return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
191 }
192
193 // Same as GetInverseElementIterator but the create iterator only return
194 // wanted type elements.
195 class SMDS_MeshNode_MyIterator:public SMDS_ElemIterator
196 {
197 private:
198   SMDS_Mesh* myMesh;
199   vtkIdType* myCells;
200   int  myNcells;
201   SMDSAbs_ElementType                                 myType;
202   int  iter;
203   vector<SMDS_MeshElement*> myFiltCells;
204
205 public:
206   SMDS_MeshNode_MyIterator(SMDS_Mesh *mesh,
207                            vtkIdType* cells,
208                            int ncells,
209                            SMDSAbs_ElementType type):
210     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
211   {
212     //MESSAGE("myNcells " << myNcells);
213     for (; iter<ncells; iter++)
214     {
215       int vtkId = myCells[iter];
216       int smdsId = myMesh->fromVtkToSmds(vtkId);
217       //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId);
218       const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
219       if (elem->GetType() == type)
220         myFiltCells.push_back((SMDS_MeshElement*)elem);
221     }
222     myNcells = myFiltCells.size();
223     //MESSAGE("myNcells " << myNcells);
224     iter = 0;
225     //MESSAGE("SMDS_MeshNode_MyIterator (filter) " << ncells << " " << myNcells);
226   }
227
228   bool more()
229   {
230     return (iter< myNcells);
231   }
232
233   const SMDS_MeshElement* next()
234   {
235     const SMDS_MeshElement* elem = myFiltCells[iter];
236     iter++;
237     return elem;
238   }
239 };
240
241 SMDS_ElemIteratorPtr SMDS_MeshNode::
242 elementsIterator(SMDSAbs_ElementType type) const
243 {
244   if(type==SMDSAbs_Node)
245     return SMDS_MeshElement::elementsIterator(SMDSAbs_Node);
246   else
247   {
248     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
249     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
250   }
251 }
252
253 int SMDS_MeshNode::NbNodes() const
254 {
255   return 1;
256 }
257
258 double* SMDS_MeshNode::getCoord() const
259 {
260   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
261 }
262
263 double SMDS_MeshNode::X() const
264 {
265   double *coord = getCoord();
266   return coord[0];
267 }
268
269 double SMDS_MeshNode::Y() const
270 {
271   double *coord = getCoord();
272   return coord[1];
273 }
274
275 double SMDS_MeshNode::Z() const
276 {
277   double *coord = getCoord();
278   return coord[2];
279 }
280
281 //================================================================================
282 /*!
283  * \brief thread safe getting coords
284  */
285 void SMDS_MeshNode::GetXYZ(double xyz[3]) const
286 {
287   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID,xyz);
288 }
289
290 //* resize the vtkPoints structure every SMDS_Mesh::chunkSize points
291 void SMDS_MeshNode::setXYZ(double x, double y, double z)
292 {
293   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
294   vtkPoints *points = mesh->getGrid()->GetPoints();
295   points->InsertPoint(myVtkID, x, y, z);
296   mesh->adjustBoundingBox(x, y, z);
297   mesh->setMyModified();
298 }
299
300 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
301 {
302   return SMDSAbs_Node;
303 }
304
305 vtkIdType SMDS_MeshNode::GetVtkType() const
306 {
307   return VTK_VERTEX;
308 }
309
310 //=======================================================================
311 //function : AddInverseElement
312 //purpose  :
313 //=======================================================================
314 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
315 {
316   const SMDS_MeshCell *cell = dynamic_cast<const SMDS_MeshCell*> (ME);
317   assert(cell);
318   SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
319   vtkCellLinks *Links = grid->GetLinks();
320   Links->ResizeCellList(myVtkID, 1);
321   Links->AddCellReference(cell->getVtkId(), myVtkID);
322 }
323
324 //=======================================================================
325 //function : ClearInverseElements
326 //purpose  :
327 //=======================================================================
328 void SMDS_MeshNode::ClearInverseElements()
329 {
330   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
331 }
332
333 bool SMDS_MeshNode::emptyInverseElements()
334 {
335   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
336   return (l.ncells == 0);
337 }
338
339 //================================================================================
340 /*!
341  * \brief Count inverse elements of given type
342  */
343 //================================================================================
344
345 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
346 {
347   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
348
349   if ( type == SMDSAbs_All )
350     return l.ncells;
351
352   int nb = 0;
353   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
354   for (int i=0; i<l.ncells; i++)
355   {
356     const SMDS_MeshElement* elem = mesh->FindElement(mesh->fromVtkToSmds(l.cells[i]));
357     if (elem->GetType() == type)
358       nb++;
359   }
360   return nb;
361 }
362
363 ///////////////////////////////////////////////////////////////////////////////
364 /// To be used with STL set
365 ///////////////////////////////////////////////////////////////////////////////
366 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
367 {
368   return e1.getVtkId()<e2.getVtkId();
369   /*if(e1.myX<e2.myX) return true;
370     else if(e1.myX==e2.myX)
371     {
372     if(e1.myY<e2.myY) return true;
373     else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);
374     else return false;
375     }
376     else return false;*/
377 }
378