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