Salome HOME
Update copyright
[modules/smesh.git] / src / SMDS / SMDS_MeshNode.cxx
1 // Copyright (C) 2007-2011  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 &&
117        myPosition != SMDS_SpacePosition::originSpacePosition() &&
118        myPosition != aPos )
119     delete myPosition;
120   myPosition = aPos;
121 }
122
123 //=======================================================================
124 //function : GetPosition
125 //purpose  : 
126 //=======================================================================
127
128 const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
129 {
130         return myPosition;
131 }
132
133 //=======================================================================
134 /*!
135  * \brief Iterator on list of elements
136  */
137 //=======================================================================
138
139 class SMDS_MeshNode_MyInvIterator: public SMDS_ElemIterator
140 {
141 private:
142   SMDS_Mesh* myMesh;
143   vtkIdType* myCells;
144   int myNcells;
145   SMDSAbs_ElementType myType;
146   int iter;
147   vector<vtkIdType> cellList;
148
149 public:
150   SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
151     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
152   {
153     //MESSAGE("SMDS_MeshNode_MyInvIterator : ncells " << myNcells);
154     cellList.clear();
155     if (type == SMDSAbs_All)
156       for (int i = 0; i < ncells; i++)
157         cellList.push_back(cells[i]);
158     else for (int i = 0; i < ncells; i++)
159       {
160         int vtkId = cells[i];
161         int smdsId = myMesh->fromVtkToSmds(vtkId);
162         const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
163         if (elem->GetType() == type)
164           {
165             //MESSAGE("Add element vtkId " << vtkId << " " << elem->GetType())
166             cellList.push_back(vtkId);
167           }
168       }
169     myCells = &cellList[0];
170     myNcells = cellList.size();
171     //MESSAGE("myNcells="<<myNcells);
172   }
173
174   bool more()
175   {
176     //MESSAGE("iter " << iter << " ncells " << myNcells);
177     return (iter < myNcells);
178   }
179
180   const SMDS_MeshElement* next()
181   {
182     int vtkId = myCells[iter];
183     int smdsId = myMesh->fromVtkToSmds(vtkId);
184     const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
185     if (!elem)
186       {
187         MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
188         throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
189       }
190     //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId << " " << elem->GetType());
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(myVtkID);
200     //MESSAGE("myID " << myID << " 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         //MESSAGE("myNcells " << myNcells);
224        for (; iter<ncells; iter++)
225         {
226            int vtkId = myCells[iter];
227            int smdsId = myMesh->fromVtkToSmds(vtkId);
228            //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId);
229            const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
230            if (elem->GetType() == type)
231                myFiltCells.push_back((SMDS_MeshElement*)elem);
232         }
233         myNcells = myFiltCells.size();
234         //MESSAGE("myNcells " << myNcells);
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(myVtkID);
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 double* SMDS_MeshNode::getCoord() const
270 {
271   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
272 }
273
274 double SMDS_MeshNode::X() const
275 {
276   double *coord = getCoord();
277   return coord[0];
278 }
279
280 double SMDS_MeshNode::Y() const
281 {
282   double *coord = getCoord();
283   return coord[1];
284 }
285
286 double SMDS_MeshNode::Z() const
287 {
288   double *coord = getCoord();
289   return coord[2];
290 }
291
292 //* resize the vtkPoints structure every SMDS_Mesh::chunkSize points
293 void SMDS_MeshNode::setXYZ(double x, double y, double z)
294 {
295   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
296   vtkPoints *points = mesh->getGrid()->GetPoints();
297   points->InsertPoint(myVtkID, x, y, z);
298   mesh->adjustBoundingBox(x, y, z);
299   mesh->setMyModified();
300 }
301
302 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
303 {
304         return SMDSAbs_Node;
305 }
306
307 vtkIdType SMDS_MeshNode::GetVtkType() const
308 {
309   return VTK_VERTEX;
310 }
311
312 //=======================================================================
313 //function : AddInverseElement
314 //purpose  :
315 //=======================================================================
316 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
317 {
318   const SMDS_MeshCell *cell = dynamic_cast<const SMDS_MeshCell*> (ME);
319   assert(cell);
320   SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
321   vtkCellLinks *Links = grid->GetCellLinks();
322   Links->ResizeCellList(myVtkID, 1);
323   Links->AddCellReference(cell->getVtkId(), myVtkID);
324 }
325
326 //=======================================================================
327 //function : ClearInverseElements
328 //purpose  :
329 //=======================================================================
330 void SMDS_MeshNode::ClearInverseElements()
331 {
332   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
333 }
334
335 bool SMDS_MeshNode::emptyInverseElements()
336 {
337   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
338   return (l.ncells == 0);
339 }
340
341 //================================================================================
342 /*!
343  * \brief Count inverse elements of given type
344  */
345 //================================================================================
346
347 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
348 {
349   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
350
351   if ( type == SMDSAbs_All )
352     return l.ncells;
353
354   int nb = 0;
355   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
356   for (int i=0; i<l.ncells; i++)
357         {
358            const SMDS_MeshElement* elem = mesh->FindElement(mesh->fromVtkToSmds(l.cells[i]));
359            if (elem->GetType() == type)
360                 nb++;
361         }
362   return nb;
363 }
364
365 ///////////////////////////////////////////////////////////////////////////////
366 /// To be used with STL set
367 ///////////////////////////////////////////////////////////////////////////////
368 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
369 {
370         return e1.getVtkId()<e2.getVtkId();
371         /*if(e1.myX<e2.myX) return true;
372         else if(e1.myX==e2.myX)
373         {
374                 if(e1.myY<e2.myY) return true;
375                 else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);
376                 else return false;
377         }
378         else return false;*/
379 }
380