Salome HOME
PR: merged from V5_1_4rc1
[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 #define protected public
36 #include <vtkCellLinks.h>
37 #define protected protected
38
39 #include "utilities.h"
40 #include "Utils_SALOME_Exception.hxx"
41 #include <cassert>
42
43 using namespace std;
44
45 int SMDS_MeshNode::nbNodes =0;
46
47 //=======================================================================
48 //function : SMDS_MeshNode
49 //purpose  : 
50 //=======================================================================
51 SMDS_MeshNode::SMDS_MeshNode() :
52   SMDS_MeshElement(-1, -1, -1),
53   myPosition(SMDS_SpacePosition::originSpacePosition())
54 {
55   nbNodes++;
56 }
57
58 SMDS_MeshNode::SMDS_MeshNode(int id, int meshId, int shapeId, double x, double y, double z):
59   SMDS_MeshElement(id, meshId, shapeId),
60   myPosition(SMDS_SpacePosition::originSpacePosition())
61 {
62   nbNodes++;
63   init(id, meshId, shapeId, x, y ,z);
64 }
65
66 void SMDS_MeshNode::init(int id, int meshId, int shapeId, double x, double y, double z)
67 {
68   myID = id;
69   myMeshId = meshId;
70   myShapeId = shapeId;
71   myIdInShape = -1;
72   //MESSAGE("Node " << myID << " (" << x << ", " << y << ", " << z << ")");
73   SMDS_Mesh* mesh = SMDS_Mesh::_meshList[myMeshId];
74   vtkUnstructuredGrid * grid = mesh->getGrid();
75   vtkPoints *points = grid->GetPoints();
76   points->InsertPoint(myID, x, y, z);
77   vtkCellLinks *cellLinks = grid->GetCellLinks();
78   if (myID >=cellLinks->Size)
79           cellLinks->Resize(myID+SMDS_Mesh::chunkSize);
80 }
81
82 SMDS_MeshNode::~SMDS_MeshNode()
83 {
84   nbNodes--;
85 }
86
87 //=======================================================================
88 //function : RemoveInverseElement
89 //purpose  : 
90 //=======================================================================
91
92 void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * parent)
93 {
94     //MESSAGE("RemoveInverseElement " << myID << " " << parent->GetID());
95     const SMDS_MeshCell* cell = dynamic_cast<const SMDS_MeshCell*>(parent);
96     MYASSERT(cell);
97     SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myID, cell->getVtkId());
98 }
99
100 //=======================================================================
101 //function : Print
102 //purpose  : 
103 //=======================================================================
104
105 void SMDS_MeshNode::Print(ostream & OS) const
106 {
107         OS << "Node <" << myID << "> : X = " << X() << " Y = "
108                 << Y() << " Z = " << Z() << endl;
109 }
110
111 //=======================================================================
112 //function : SetPosition
113 //purpose  : 
114 //=======================================================================
115
116 void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
117 {
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
146  public:
147   SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh,
148                               vtkIdType* cells,
149                               int ncells,
150                               SMDSAbs_ElementType type):
151     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
152   {
153      //MESSAGE("SMDS_MeshNode_MyInvIterator : ncells " << myNcells);
154   }
155
156   bool more()
157   {
158       return (iter< myNcells);
159   }
160
161   const SMDS_MeshElement* next()
162   {
163       int vtkId = myCells[iter];
164       int smdsId = myMesh->fromVtkToSmds(vtkId);
165       const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
166       if (!elem)
167       {
168           assert(0);
169           throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
170       }
171       //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId << " " << (elem!=0));
172       iter++;
173       return elem;
174   }
175 };
176
177 SMDS_ElemIteratorPtr SMDS_MeshNode::
178         GetInverseElementIterator(SMDSAbs_ElementType type) const
179 {
180     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myID);
181     //MESSAGE("myID " << myID << " ncells " << l.ncells);
182     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
183 }
184
185 // Same as GetInverseElementIterator but the create iterator only return
186 // wanted type elements.
187 class SMDS_MeshNode_MyIterator:public SMDS_ElemIterator
188 {
189 private:
190   SMDS_Mesh* myMesh;
191   vtkIdType* myCells;
192   int  myNcells;
193   SMDSAbs_ElementType                                 myType;
194   int  iter;
195   vector<SMDS_MeshElement*> myFiltCells;
196
197  public:
198   SMDS_MeshNode_MyIterator(SMDS_Mesh *mesh,
199                            vtkIdType* cells,
200                            int ncells,
201                            SMDSAbs_ElementType type):
202     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
203   {
204         //MESSAGE("myNcells " << myNcells);
205        for (; iter<ncells; iter++)
206         {
207            int vtkId = myCells[iter];
208            int smdsId = myMesh->fromVtkToSmds(vtkId);
209            //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId);
210            const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
211            if (elem->GetType() == type)
212                myFiltCells.push_back((SMDS_MeshElement*)elem);
213         }
214         myNcells = myFiltCells.size();
215         //MESSAGE("myNcells " << myNcells);
216        iter = 0;
217         //MESSAGE("SMDS_MeshNode_MyIterator (filter) " << ncells << " " << myNcells);
218   }
219
220   bool more()
221   {
222       return (iter< myNcells);
223   }
224
225   const SMDS_MeshElement* next()
226   {
227       const SMDS_MeshElement* elem = myFiltCells[iter];
228       iter++;
229       return elem;
230   }
231 };
232
233 SMDS_ElemIteratorPtr SMDS_MeshNode::
234         elementsIterator(SMDSAbs_ElementType type) const
235 {
236   if(type==SMDSAbs_Node)
237     return SMDS_MeshElement::elementsIterator(SMDSAbs_Node); 
238   else
239   {
240     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myID);
241     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
242   }
243 }
244
245 int SMDS_MeshNode::NbNodes() const
246 {
247         return 1;
248 }
249
250
251 double* SMDS_MeshNode::getCoord() const
252 {
253   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myID);
254 }
255
256 double SMDS_MeshNode::X() const
257 {
258   double *coord = getCoord();
259   return coord[0];
260 }
261
262 double SMDS_MeshNode::Y() const
263 {
264   double *coord = getCoord();
265   return coord[1];
266 }
267
268 double SMDS_MeshNode::Z() const
269 {
270   double *coord = getCoord();
271   return coord[2];
272 }
273
274 //* resize the vtkPoints structure every SMDS_Mesh::chunkSize points
275 void SMDS_MeshNode::setXYZ(double x, double y, double z)
276 {
277   vtkPoints *points = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoints();
278   points->InsertPoint(myID, x, y, z);
279 }
280
281 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
282 {
283         return SMDSAbs_Node;
284 }
285
286 vtkIdType SMDS_MeshNode::GetVtkType() const
287 {
288   return VTK_VERTEX;
289 }
290
291 //=======================================================================
292 //function : AddInverseElement
293 //purpose  :
294 //=======================================================================
295 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
296 {
297     const SMDS_MeshCell *cell = dynamic_cast<const SMDS_MeshCell*>(ME);
298     assert(cell);
299     SMDS_Mesh::_meshList[myMeshId]->getGrid()->AddReferenceToCell(myID, cell->getVtkId());
300 }
301
302 //=======================================================================
303 //function : ClearInverseElements
304 //purpose  :
305 //=======================================================================
306 void SMDS_MeshNode::ClearInverseElements()
307 {
308   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myID, 0);
309 }
310
311 bool SMDS_MeshNode::emptyInverseElements()
312 {
313   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myID);
314   return (l.ncells == 0);
315 }
316
317 //================================================================================
318 /*!
319  * \brief Count inverse elements of given type
320  */
321 //================================================================================
322
323 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
324 {
325   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myID);
326
327   if ( type == SMDSAbs_All )
328     return l.ncells;
329
330   int nb = 0;
331   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
332   for (int i=0; i<l.ncells; i++)
333         {
334            const SMDS_MeshElement* elem = mesh->FindNode(l.cells[i]);
335            if (elem->GetType() == type)
336                 nb++;
337         }
338   return nb;
339 }
340
341 ///////////////////////////////////////////////////////////////////////////////
342 /// To be used with STL set
343 ///////////////////////////////////////////////////////////////////////////////
344 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
345 {
346         return e1.GetID()<e2.GetID();
347         /*if(e1.myX<e2.myX) return true;
348         else if(e1.myX==e2.myX)
349         {
350                 if(e1.myY<e2.myY) return true;
351                 else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);
352                 else return false;
353         }
354         else return false;*/
355 }
356