Salome HOME
Regression of doc/salome/examples/transforming_meshes_ex10.py
[modules/smesh.git] / src / SMDS / SMDS_MeshNode.cxx
1 // Copyright (C) 2007-2016  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, or (at your option) any later version.
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   SMDS_UnstructuredGrid * grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
68   vtkPoints *points = grid->GetPoints();
69   points->InsertPoint(myVtkID, x, y, z);
70   if ( grid->HasLinks() )
71     grid->GetLinks()->ResizeForPoint( myVtkID );
72 }
73
74 SMDS_MeshNode::~SMDS_MeshNode()
75 {
76   nbNodes--;
77   if ( myPosition && myPosition != SMDS_SpacePosition::originSpacePosition() )
78     delete myPosition, myPosition = 0;
79 }
80
81 //=======================================================================
82 //function : RemoveInverseElement
83 //purpose  :
84 //=======================================================================
85
86 void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * elem)
87 {
88   if ( SMDS_Mesh::_meshList[myMeshId]->getGrid()->HasLinks() )
89     SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, elem->getVtkId());
90 }
91
92 //=======================================================================
93 //function : Print
94 //purpose  :
95 //=======================================================================
96
97 void SMDS_MeshNode::Print(ostream & OS) const
98 {
99   OS << "Node <" << myID << "> : X = " << X() << " Y = "
100      << Y() << " Z = " << Z() << endl;
101 }
102
103 //=======================================================================
104 //function : SetPosition
105 //purpose  :
106 //=======================================================================
107
108 void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
109 {
110   if ( myPosition &&
111        myPosition != SMDS_SpacePosition::originSpacePosition() &&
112        myPosition != aPos )
113     delete myPosition;
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     if ( ncells )
148     {
149       cellList.reserve( ncells );
150       if (type == SMDSAbs_All)
151       {
152         cellList.assign( cells, cells + ncells );
153       }
154       else
155       {
156         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             cellList.push_back(vtkId);
164           }
165         }
166       }
167       myCells = cellList.empty() ? 0 : &cellList[0];
168       myNcells = cellList.size();
169     }
170   }
171
172   bool more()
173   {
174     return (iter < myNcells);
175   }
176
177   const SMDS_MeshElement* next()
178   {
179     int vtkId = myCells[iter];
180     int smdsId = myMesh->fromVtkToSmds(vtkId);
181     const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
182     if (!elem)
183     {
184       MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
185       throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
186     }
187     iter++;
188     return elem;
189   }
190 };
191
192 SMDS_ElemIteratorPtr SMDS_MeshNode::GetInverseElementIterator(SMDSAbs_ElementType type) const
193 {
194   if ( SMDS_Mesh::_meshList[myMeshId]->NbElements() > 0 ) // avoid building links
195   {
196     vtkCellLinks::Link& l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
197     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
198   }
199   else
200   {
201     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], 0, 0, type));
202   }
203 }
204
205 SMDS_ElemIteratorPtr SMDS_MeshNode::elementsIterator(SMDSAbs_ElementType type) const
206 {
207   if ( type == SMDSAbs_Node )
208     return SMDS_MeshElement::elementsIterator( SMDSAbs_Node );
209   else
210     return GetInverseElementIterator( type );
211 }
212
213 int SMDS_MeshNode::NbNodes() const
214 {
215   return 1;
216 }
217
218 double* SMDS_MeshNode::getCoord() const
219 {
220   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
221 }
222
223 double SMDS_MeshNode::X() const
224 {
225   double *coord = getCoord();
226   return coord[0];
227 }
228
229 double SMDS_MeshNode::Y() const
230 {
231   double *coord = getCoord();
232   return coord[1];
233 }
234
235 double SMDS_MeshNode::Z() const
236 {
237   double *coord = getCoord();
238   return coord[2];
239 }
240
241 //================================================================================
242 /*!
243  * \brief thread safe getting coords
244  */
245 //================================================================================
246
247 void SMDS_MeshNode::GetXYZ(double xyz[3]) const
248 {
249   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID,xyz);
250 }
251
252 //================================================================================
253 void SMDS_MeshNode::setXYZ(double x, double y, double z)
254 {
255   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
256   vtkPoints *points = mesh->getGrid()->GetPoints();
257   points->InsertPoint(myVtkID, x, y, z);
258   mesh->adjustBoundingBox(x, y, z);
259   mesh->setMyModified();
260 }
261
262 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
263 {
264   return SMDSAbs_Node;
265 }
266
267 vtkIdType SMDS_MeshNode::GetVtkType() const
268 {
269   return VTK_VERTEX;
270 }
271
272 //=======================================================================
273 //function : AddInverseElement
274 //purpose  :
275 //=======================================================================
276 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
277 {
278   SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
279   if ( grid->HasLinks() )
280   {
281     vtkCellLinks *Links = grid->GetLinks();
282     Links->ResizeCellList(myVtkID, 1);
283     Links->AddCellReference(ME->getVtkId(), myVtkID);
284   }
285 }
286
287 //=======================================================================
288 //function : ClearInverseElements
289 //purpose  :
290 //=======================================================================
291 void SMDS_MeshNode::ClearInverseElements()
292 {
293   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
294 }
295
296 //================================================================================
297 /*!
298  * \brief Count inverse elements of given type
299  */
300 //================================================================================
301
302 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
303 {
304   int nb = 0;
305   if ( SMDS_Mesh::_meshList[myMeshId]->NbElements() > 0 ) // avoid building links
306   {
307     vtkCellLinks::Link& l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
308
309     if ( type == SMDSAbs_All )
310       return l.ncells;
311
312     SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
313     for ( int i = 0; i < l.ncells; i++ )
314     {
315       const SMDS_MeshElement* elem = mesh->FindElement( mesh->fromVtkToSmds( l.cells[i] ));
316       nb += ( elem->GetType() == type );
317     }
318   }
319   return nb;
320 }