Salome HOME
Fix regressions
[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 * elem)
86 {
87   if ( SMDS_Mesh::_meshList[myMeshId]->getGrid()->HasLinks() )
88     SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, elem->getVtkId());
89 }
90
91 //=======================================================================
92 //function : Print
93 //purpose  :
94 //=======================================================================
95
96 void SMDS_MeshNode::Print(ostream & OS) const
97 {
98   OS << "Node <" << myID << "> : X = " << X() << " Y = "
99      << Y() << " Z = " << Z() << endl;
100 }
101
102 //=======================================================================
103 //function : SetPosition
104 //purpose  :
105 //=======================================================================
106
107 void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
108 {
109   if ( myPosition &&
110        myPosition != SMDS_SpacePosition::originSpacePosition() &&
111        myPosition != aPos )
112     delete myPosition;
113   myPosition = aPos;
114 }
115
116 //=======================================================================
117 //function : GetPosition
118 //purpose  :
119 //=======================================================================
120
121 const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
122 {
123   return myPosition;
124 }
125
126 //=======================================================================
127 /*!
128  * \brief Iterator on list of elements
129  */
130 //=======================================================================
131
132 class SMDS_MeshNode_MyInvIterator: public SMDS_ElemIterator
133 {
134 private:
135   SMDS_Mesh* myMesh;
136   vtkIdType* myCells;
137   int myNcells;
138   SMDSAbs_ElementType myType;
139   int iter;
140   vector<vtkIdType> cellList;
141
142 public:
143   SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
144     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
145   {
146     cellList.reserve( ncells );
147     if (type == SMDSAbs_All)
148       cellList.assign( cells, cells + ncells );
149     else
150       for (int i = 0; i < ncells; i++)
151       {
152         int vtkId = cells[i];
153         int smdsId = myMesh->fromVtkToSmds(vtkId);
154         const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
155         if (elem->GetType() == type)
156         {
157           cellList.push_back(vtkId);
158         }
159       }
160     myCells = cellList.empty() ? 0 : &cellList[0];
161     myNcells = cellList.size();
162   }
163
164   bool more()
165   {
166     return (iter < myNcells);
167   }
168
169   const SMDS_MeshElement* next()
170   {
171     int vtkId = myCells[iter];
172     int smdsId = myMesh->fromVtkToSmds(vtkId);
173     const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
174     if (!elem)
175     {
176       MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
177       throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
178     }
179     iter++;
180     return elem;
181   }
182 };
183
184 SMDS_ElemIteratorPtr SMDS_MeshNode::GetInverseElementIterator(SMDSAbs_ElementType type) const
185 {
186   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
187   return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
188 }
189
190 // Same as GetInverseElementIterator but the created iterator only returns
191 // wanted type elements.
192 class SMDS_MeshNode_MyIterator:public SMDS_ElemIterator
193 {
194 private:
195   SMDS_Mesh*                myMesh;
196   vtkIdType*                myCells;
197   int                       myNcells;
198   SMDSAbs_ElementType       myType;
199   int                       iter;
200   vector<SMDS_MeshElement*> myFiltCells;
201
202 public:
203   SMDS_MeshNode_MyIterator(SMDS_Mesh *mesh,
204                            vtkIdType* cells,
205                            int ncells,
206                            SMDSAbs_ElementType type):
207     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
208   {
209     for (; iter<ncells; iter++)
210     {
211       int vtkId = myCells[iter];
212       int smdsId = myMesh->fromVtkToSmds(vtkId);
213       const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
214       if (elem->GetType() == type)
215         myFiltCells.push_back((SMDS_MeshElement*)elem);
216     }
217     myNcells = myFiltCells.size();
218     iter = 0;
219   }
220
221   bool more()
222   {
223     return (iter< myNcells);
224   }
225
226   const SMDS_MeshElement* next()
227   {
228     const SMDS_MeshElement* elem = myFiltCells[iter];
229     iter++;
230     return elem;
231   }
232 };
233
234 SMDS_ElemIteratorPtr SMDS_MeshNode::
235 elementsIterator(SMDSAbs_ElementType type) const
236 {
237   if(type==SMDSAbs_Node)
238     return SMDS_MeshElement::elementsIterator(SMDSAbs_Node);
239   else
240   {
241     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
242     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
243   }
244 }
245
246 int SMDS_MeshNode::NbNodes() const
247 {
248   return 1;
249 }
250
251 double* SMDS_MeshNode::getCoord() const
252 {
253   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
254 }
255
256 double SMDS_MeshNode::X() const
257 {
258   double *coord = getCoord();
259   return coord[0];
260 }
261
262 double SMDS_MeshNode::Y() const
263 {
264   double *coord = getCoord();
265   return coord[1];
266 }
267
268 double SMDS_MeshNode::Z() const
269 {
270   double *coord = getCoord();
271   return coord[2];
272 }
273
274 //================================================================================
275 /*!
276  * \brief thread safe getting coords
277  */
278 void SMDS_MeshNode::GetXYZ(double xyz[3]) const
279 {
280   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID,xyz);
281 }
282
283 //* resize the vtkPoints structure every SMDS_Mesh::chunkSize points
284 void SMDS_MeshNode::setXYZ(double x, double y, double z)
285 {
286   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
287   vtkPoints *points = mesh->getGrid()->GetPoints();
288   points->InsertPoint(myVtkID, x, y, z);
289   mesh->adjustBoundingBox(x, y, z);
290   mesh->setMyModified();
291 }
292
293 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
294 {
295   return SMDSAbs_Node;
296 }
297
298 vtkIdType SMDS_MeshNode::GetVtkType() const
299 {
300   return VTK_VERTEX;
301 }
302
303 //=======================================================================
304 //function : AddInverseElement
305 //purpose  :
306 //=======================================================================
307 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
308 {
309   SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
310   vtkCellLinks *Links = grid->GetLinks();
311   Links->ResizeCellList(myVtkID, 1);
312   Links->AddCellReference(ME->getVtkId(), myVtkID);
313 }
314
315 //=======================================================================
316 //function : ClearInverseElements
317 //purpose  :
318 //=======================================================================
319 void SMDS_MeshNode::ClearInverseElements()
320 {
321   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
322 }
323
324 bool SMDS_MeshNode::emptyInverseElements()
325 {
326   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
327   return (l.ncells == 0);
328 }
329
330 //================================================================================
331 /*!
332  * \brief Count inverse elements of given type
333  */
334 //================================================================================
335
336 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
337 {
338   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
339
340   if ( type == SMDSAbs_All )
341     return l.ncells;
342
343   int nb = 0;
344   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
345   for (int i=0; i<l.ncells; i++)
346   {
347     const SMDS_MeshElement* elem = mesh->FindElement(mesh->fromVtkToSmds(l.cells[i]));
348     if (elem->GetType() == type)
349       nb++;
350   }
351   return nb;
352 }
353
354 ///////////////////////////////////////////////////////////////////////////////
355 /// To be used with STL set
356 ///////////////////////////////////////////////////////////////////////////////
357 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
358 {
359   return e1.getVtkId()<e2.getVtkId();
360   /*if(e1.myX<e2.myX) return true;
361     else if(e1.myX==e2.myX)
362     {
363     if(e1.myY<e2.myY) return true;
364     else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);
365     else return false;
366     }
367     else return false;*/
368 }
369