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