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