1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
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.
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.
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
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // SMESH SMDS : implementaion of Salome mesh data structure
26 #pragma warning(disable:4786)
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_SpacePosition.hxx"
31 #include "SMDS_IteratorOfElements.hxx"
32 #include "SMDS_Mesh.hxx"
33 #include <vtkUnstructuredGrid.h>
35 #include "utilities.h"
36 #include "Utils_SALOME_Exception.hxx"
41 int SMDS_MeshNode::nbNodes =0;
43 //=======================================================================
44 //function : SMDS_MeshNode
46 //=======================================================================
47 SMDS_MeshNode::SMDS_MeshNode() :
48 SMDS_MeshElement(-1, -1, 0),
49 myPosition(SMDS_SpacePosition::originSpacePosition())
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())
59 init(id, meshId, shapeId, x, y ,z);
62 void SMDS_MeshNode::init(int id, int meshId, int shapeId, double x, double y, double z)
64 SMDS_MeshElement::init(id, meshId, shapeId);
67 SMDS_UnstructuredGrid * grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
68 vtkPoints *points = grid->GetPoints();
69 points->InsertPoint(myVtkID, x, y, z);
70 grid->GetLinks()->ResizeForPoint( myVtkID );
73 SMDS_MeshNode::~SMDS_MeshNode()
76 if ( myPosition && myPosition != SMDS_SpacePosition::originSpacePosition() )
77 delete myPosition, myPosition = 0;
80 //=======================================================================
81 //function : RemoveInverseElement
83 //=======================================================================
85 void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * elem)
87 if ( SMDS_Mesh::_meshList[myMeshId]->getGrid()->HasLinks() )
88 SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, elem->getVtkId());
91 //=======================================================================
94 //=======================================================================
96 void SMDS_MeshNode::Print(ostream & OS) const
98 OS << "Node <" << myID << "> : X = " << X() << " Y = "
99 << Y() << " Z = " << Z() << endl;
102 //=======================================================================
103 //function : SetPosition
105 //=======================================================================
107 void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
110 myPosition != SMDS_SpacePosition::originSpacePosition() &&
116 //=======================================================================
117 //function : GetPosition
119 //=======================================================================
121 const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
126 //=======================================================================
128 * \brief Iterator on list of elements
130 //=======================================================================
132 class SMDS_MeshNode_MyInvIterator: public SMDS_ElemIterator
138 SMDSAbs_ElementType myType;
140 vector<vtkIdType> cellList;
143 SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
144 myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
146 cellList.reserve( ncells );
147 if (type == SMDSAbs_All)
148 cellList.assign( cells, cells + ncells );
150 for (int i = 0; i < ncells; i++)
152 int vtkId = cells[i];
153 int smdsId = myMesh->fromVtkToSmds(vtkId);
154 const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
155 if (elem->GetType() == type)
157 cellList.push_back(vtkId);
160 myCells = cellList.empty() ? 0 : &cellList[0];
161 myNcells = cellList.size();
166 return (iter < myNcells);
169 const SMDS_MeshElement* next()
171 int vtkId = myCells[iter];
172 int smdsId = myMesh->fromVtkToSmds(vtkId);
173 const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
176 MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
177 throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
184 SMDS_ElemIteratorPtr SMDS_MeshNode::GetInverseElementIterator(SMDSAbs_ElementType type) const
186 vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
187 return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
190 // Same as GetInverseElementIterator but the created iterator only returns
191 // wanted type elements.
192 class SMDS_MeshNode_MyIterator:public SMDS_ElemIterator
198 SMDSAbs_ElementType myType;
200 vector<SMDS_MeshElement*> myFiltCells;
203 SMDS_MeshNode_MyIterator(SMDS_Mesh *mesh,
206 SMDSAbs_ElementType type):
207 myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
209 for (; iter<ncells; iter++)
211 int vtkId = myCells[iter];
212 int smdsId = myMesh->fromVtkToSmds(vtkId);
213 const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
214 if (elem->GetType() == type)
215 myFiltCells.push_back((SMDS_MeshElement*)elem);
217 myNcells = myFiltCells.size();
223 return (iter< myNcells);
226 const SMDS_MeshElement* next()
228 const SMDS_MeshElement* elem = myFiltCells[iter];
234 SMDS_ElemIteratorPtr SMDS_MeshNode::
235 elementsIterator(SMDSAbs_ElementType type) const
237 if(type==SMDSAbs_Node)
238 return SMDS_MeshElement::elementsIterator(SMDSAbs_Node);
241 vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
242 return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
246 int SMDS_MeshNode::NbNodes() const
251 double* SMDS_MeshNode::getCoord() const
253 return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
256 double SMDS_MeshNode::X() const
258 double *coord = getCoord();
262 double SMDS_MeshNode::Y() const
264 double *coord = getCoord();
268 double SMDS_MeshNode::Z() const
270 double *coord = getCoord();
274 //================================================================================
276 * \brief thread safe getting coords
278 void SMDS_MeshNode::GetXYZ(double xyz[3]) const
280 return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID,xyz);
283 //* resize the vtkPoints structure every SMDS_Mesh::chunkSize points
284 void SMDS_MeshNode::setXYZ(double x, double y, double z)
286 SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
287 vtkPoints *points = mesh->getGrid()->GetPoints();
288 points->InsertPoint(myVtkID, x, y, z);
289 mesh->adjustBoundingBox(x, y, z);
290 mesh->setMyModified();
293 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
298 vtkIdType SMDS_MeshNode::GetVtkType() const
303 //=======================================================================
304 //function : AddInverseElement
306 //=======================================================================
307 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
309 SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
310 vtkCellLinks *Links = grid->GetLinks();
311 Links->ResizeCellList(myVtkID, 1);
312 Links->AddCellReference(ME->getVtkId(), myVtkID);
315 //=======================================================================
316 //function : ClearInverseElements
318 //=======================================================================
319 void SMDS_MeshNode::ClearInverseElements()
321 SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
324 bool SMDS_MeshNode::emptyInverseElements()
326 vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
327 return (l.ncells == 0);
330 //================================================================================
332 * \brief Count inverse elements of given type
334 //================================================================================
336 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
338 vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetLinks()->GetLink(myVtkID);
340 if ( type == SMDSAbs_All )
344 SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
345 for (int i=0; i<l.ncells; i++)
347 const SMDS_MeshElement* elem = mesh->FindElement(mesh->fromVtkToSmds(l.cells[i]));
348 if (elem->GetType() == type)
354 ///////////////////////////////////////////////////////////////////////////////
355 /// To be used with STL set
356 ///////////////////////////////////////////////////////////////////////////////
357 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
359 return e1.getVtkId()<e2.getVtkId();
360 /*if(e1.myX<e2.myX) return true;
361 else if(e1.myX==e2.myX)
363 if(e1.myY<e2.myY) return true;
364 else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);