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