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