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