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