]> SALOME platform Git repositories - modules/smesh.git/blob - src/SMDS/SMDS_MeshNode.cxx
Salome HOME
delete node positions at node removal
[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   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   if (myVtkID >= cellLinks->GetLinksSize())
75           cellLinks->ResizeL(myVtkID+SMDS_Mesh::chunkSize);
76 }
77
78 SMDS_MeshNode::~SMDS_MeshNode()
79 {
80   nbNodes--;
81   if ( myPosition && myPosition != SMDS_SpacePosition::originSpacePosition() )
82     delete myPosition, myPosition = 0;
83 }
84
85 //=======================================================================
86 //function : RemoveInverseElement
87 //purpose  : 
88 //=======================================================================
89
90 void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * parent)
91 {
92     //MESSAGE("RemoveInverseElement " << myID << " " << parent->GetID());
93     const SMDS_MeshCell* cell = dynamic_cast<const SMDS_MeshCell*>(parent);
94     MYASSERT(cell);
95     SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, cell->getVtkId());
96 }
97
98 //=======================================================================
99 //function : Print
100 //purpose  : 
101 //=======================================================================
102
103 void SMDS_MeshNode::Print(ostream & OS) const
104 {
105         OS << "Node <" << myID << "> : X = " << X() << " Y = "
106                 << Y() << " Z = " << Z() << endl;
107 }
108
109 //=======================================================================
110 //function : SetPosition
111 //purpose  : 
112 //=======================================================================
113
114 void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
115 {
116   if ( myPosition && myPosition != aPos )
117     delete myPosition;
118   myPosition = aPos;
119 }
120
121 //=======================================================================
122 //function : GetPosition
123 //purpose  : 
124 //=======================================================================
125
126 const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
127 {
128         return myPosition;
129 }
130
131 //=======================================================================
132 /*!
133  * \brief Iterator on list of elements
134  */
135 //=======================================================================
136
137 class SMDS_MeshNode_MyInvIterator: public SMDS_ElemIterator
138 {
139 private:
140   SMDS_Mesh* myMesh;
141   vtkIdType* myCells;
142   int myNcells;
143   SMDSAbs_ElementType myType;
144   int iter;
145   vector<vtkIdType> cellList;
146
147 public:
148   SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
149     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
150   {
151     //MESSAGE("SMDS_MeshNode_MyInvIterator : ncells " << myNcells);
152     cellList.clear();
153     if (type == SMDSAbs_All)
154       for (int i = 0; i < ncells; i++)
155         cellList.push_back(cells[i]);
156     else 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             //MESSAGE("Add element vtkId " << vtkId << " " << elem->GetType())
164             cellList.push_back(vtkId);
165           }
166       }
167     myCells = &cellList[0];
168     myNcells = cellList.size();
169     //MESSAGE("myNcells="<<myNcells);
170   }
171
172   bool more()
173   {
174     //MESSAGE("iter " << iter << " ncells " << myNcells);
175     return (iter < myNcells);
176   }
177
178   const SMDS_MeshElement* next()
179   {
180     int vtkId = myCells[iter];
181     int smdsId = myMesh->fromVtkToSmds(vtkId);
182     const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
183     if (!elem)
184       {
185         MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
186         throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
187       }
188     //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId << " " << elem->GetType());
189     iter++;
190     return elem;
191   }
192 };
193
194 SMDS_ElemIteratorPtr SMDS_MeshNode::
195         GetInverseElementIterator(SMDSAbs_ElementType type) const
196 {
197     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
198     //MESSAGE("myID " << myID << " ncells " << l.ncells);
199     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
200 }
201
202 // Same as GetInverseElementIterator but the create iterator only return
203 // wanted type elements.
204 class SMDS_MeshNode_MyIterator:public SMDS_ElemIterator
205 {
206 private:
207   SMDS_Mesh* myMesh;
208   vtkIdType* myCells;
209   int  myNcells;
210   SMDSAbs_ElementType                                 myType;
211   int  iter;
212   vector<SMDS_MeshElement*> myFiltCells;
213
214  public:
215   SMDS_MeshNode_MyIterator(SMDS_Mesh *mesh,
216                            vtkIdType* cells,
217                            int ncells,
218                            SMDSAbs_ElementType type):
219     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
220   {
221         //MESSAGE("myNcells " << myNcells);
222        for (; iter<ncells; iter++)
223         {
224            int vtkId = myCells[iter];
225            int smdsId = myMesh->fromVtkToSmds(vtkId);
226            //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId);
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         //MESSAGE("myNcells " << myNcells);
233        iter = 0;
234         //MESSAGE("SMDS_MeshNode_MyIterator (filter) " << ncells << " " << myNcells);
235   }
236
237   bool more()
238   {
239       return (iter< myNcells);
240   }
241
242   const SMDS_MeshElement* next()
243   {
244       const SMDS_MeshElement* elem = myFiltCells[iter];
245       iter++;
246       return elem;
247   }
248 };
249
250 SMDS_ElemIteratorPtr SMDS_MeshNode::
251         elementsIterator(SMDSAbs_ElementType type) const
252 {
253   if(type==SMDSAbs_Node)
254     return SMDS_MeshElement::elementsIterator(SMDSAbs_Node); 
255   else
256   {
257     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
258     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
259   }
260 }
261
262 int SMDS_MeshNode::NbNodes() const
263 {
264         return 1;
265 }
266
267 double* SMDS_MeshNode::getCoord() const
268 {
269   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
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   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
294   vtkPoints *points = mesh->getGrid()->GetPoints();
295   points->InsertPoint(myVtkID, x, y, z);
296   mesh->adjustBoundingBox(x, y, z);
297   mesh->setMyModified();
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_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
319   vtkCellLinks *Links = grid->GetCellLinks();
320   Links->ResizeCellList(myVtkID, 1);
321   Links->AddCellReference(cell->getVtkId(), myVtkID);
322 }
323
324 //=======================================================================
325 //function : ClearInverseElements
326 //purpose  :
327 //=======================================================================
328 void SMDS_MeshNode::ClearInverseElements()
329 {
330   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
331 }
332
333 bool SMDS_MeshNode::emptyInverseElements()
334 {
335   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
336   return (l.ncells == 0);
337 }
338
339 //================================================================================
340 /*!
341  * \brief Count inverse elements of given type
342  */
343 //================================================================================
344
345 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
346 {
347   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
348
349   if ( type == SMDSAbs_All )
350     return l.ncells;
351
352   int nb = 0;
353   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
354   for (int i=0; i<l.ncells; i++)
355         {
356            const SMDS_MeshElement* elem = mesh->FindElement(mesh->fromVtkToSmds(l.cells[i]));
357            if (elem->GetType() == type)
358                 nb++;
359         }
360   return nb;
361 }
362
363 ///////////////////////////////////////////////////////////////////////////////
364 /// To be used with STL set
365 ///////////////////////////////////////////////////////////////////////////////
366 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
367 {
368         return e1.getVtkId()<e2.getVtkId();
369         /*if(e1.myX<e2.myX) return true;
370         else if(e1.myX==e2.myX)
371         {
372                 if(e1.myY<e2.myY) return true;
373                 else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);
374                 else return false;
375         }
376         else return false;*/
377 }
378