Salome HOME
23080: [CEA 1497] Do not merge a middle node in quadratic with the extreme nodes...
[modules/smesh.git] / src / SMDS / SMDS_MeshNode.cxx
1 // Copyright (C) 2007-2015  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   //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     cellList.reserve( ncells );
153     if (type == SMDSAbs_All)
154       cellList.assign( cells, cells + ncells );
155     else
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     myCells = cellList.empty() ? 0 : &cellList[0];
167     myNcells = cellList.size();
168   }
169
170   bool more()
171   {
172     return (iter < myNcells);
173   }
174
175   const SMDS_MeshElement* next()
176   {
177     int vtkId = myCells[iter];
178     int smdsId = myMesh->fromVtkToSmds(vtkId);
179     const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
180     if (!elem)
181       {
182         MESSAGE("SMDS_MeshNode_MyInvIterator problem Null element");
183         throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
184       }
185     //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId << " " << elem->GetType());
186     iter++;
187     return elem;
188   }
189 };
190
191 SMDS_ElemIteratorPtr SMDS_MeshNode::
192         GetInverseElementIterator(SMDSAbs_ElementType type) const
193 {
194     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
195     //MESSAGE("myID " << myID << " ncells " << l.ncells);
196     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyInvIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
197 }
198
199 // Same as GetInverseElementIterator but the create iterator only return
200 // wanted type elements.
201 class SMDS_MeshNode_MyIterator:public SMDS_ElemIterator
202 {
203 private:
204   SMDS_Mesh* myMesh;
205   vtkIdType* myCells;
206   int  myNcells;
207   SMDSAbs_ElementType                                 myType;
208   int  iter;
209   vector<SMDS_MeshElement*> myFiltCells;
210
211  public:
212   SMDS_MeshNode_MyIterator(SMDS_Mesh *mesh,
213                            vtkIdType* cells,
214                            int ncells,
215                            SMDSAbs_ElementType type):
216     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
217   {
218         //MESSAGE("myNcells " << myNcells);
219        for (; iter<ncells; iter++)
220         {
221            int vtkId = myCells[iter];
222            int smdsId = myMesh->fromVtkToSmds(vtkId);
223            //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId);
224            const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
225            if (elem->GetType() == type)
226                myFiltCells.push_back((SMDS_MeshElement*)elem);
227         }
228         myNcells = myFiltCells.size();
229         //MESSAGE("myNcells " << myNcells);
230        iter = 0;
231         //MESSAGE("SMDS_MeshNode_MyIterator (filter) " << ncells << " " << myNcells);
232   }
233
234   bool more()
235   {
236       return (iter< myNcells);
237   }
238
239   const SMDS_MeshElement* next()
240   {
241       const SMDS_MeshElement* elem = myFiltCells[iter];
242       iter++;
243       return elem;
244   }
245 };
246
247 SMDS_ElemIteratorPtr SMDS_MeshNode::
248         elementsIterator(SMDSAbs_ElementType type) const
249 {
250   if(type==SMDSAbs_Node)
251     return SMDS_MeshElement::elementsIterator(SMDSAbs_Node); 
252   else
253   {
254     vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
255     return SMDS_ElemIteratorPtr(new SMDS_MeshNode_MyIterator(SMDS_Mesh::_meshList[myMeshId], l.cells, l.ncells, type));
256   }
257 }
258
259 int SMDS_MeshNode::NbNodes() const
260 {
261         return 1;
262 }
263
264 double* SMDS_MeshNode::getCoord() const
265 {
266   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID);
267 }
268
269 double SMDS_MeshNode::X() const
270 {
271   double *coord = getCoord();
272   return coord[0];
273 }
274
275 double SMDS_MeshNode::Y() const
276 {
277   double *coord = getCoord();
278   return coord[1];
279 }
280
281 double SMDS_MeshNode::Z() const
282 {
283   double *coord = getCoord();
284   return coord[2];
285 }
286
287 //================================================================================
288 /*!
289  * \brief thread safe getting coords
290  */
291 void SMDS_MeshNode::GetXYZ(double xyz[3]) const
292 {
293   return SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetPoint(myVtkID,xyz);
294 }
295
296 //* resize the vtkPoints structure every SMDS_Mesh::chunkSize points
297 void SMDS_MeshNode::setXYZ(double x, double y, double z)
298 {
299   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
300   vtkPoints *points = mesh->getGrid()->GetPoints();
301   points->InsertPoint(myVtkID, x, y, z);
302   mesh->adjustBoundingBox(x, y, z);
303   mesh->setMyModified();
304 }
305
306 SMDSAbs_ElementType SMDS_MeshNode::GetType() const
307 {
308         return SMDSAbs_Node;
309 }
310
311 vtkIdType SMDS_MeshNode::GetVtkType() const
312 {
313   return VTK_VERTEX;
314 }
315
316 //=======================================================================
317 //function : AddInverseElement
318 //purpose  :
319 //=======================================================================
320 void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
321 {
322   const SMDS_MeshCell *cell = dynamic_cast<const SMDS_MeshCell*> (ME);
323   assert(cell);
324   SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
325   vtkCellLinks *Links = grid->GetCellLinks();
326   Links->ResizeCellList(myVtkID, 1);
327   Links->AddCellReference(cell->getVtkId(), myVtkID);
328 }
329
330 //=======================================================================
331 //function : ClearInverseElements
332 //purpose  :
333 //=======================================================================
334 void SMDS_MeshNode::ClearInverseElements()
335 {
336   SMDS_Mesh::_meshList[myMeshId]->getGrid()->ResizeCellList(myVtkID, 0);
337 }
338
339 bool SMDS_MeshNode::emptyInverseElements()
340 {
341   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
342   return (l.ncells == 0);
343 }
344
345 //================================================================================
346 /*!
347  * \brief Count inverse elements of given type
348  */
349 //================================================================================
350
351 int SMDS_MeshNode::NbInverseElements(SMDSAbs_ElementType type) const
352 {
353   vtkCellLinks::Link l = SMDS_Mesh::_meshList[myMeshId]->getGrid()->GetCellLinks()->GetLink(myVtkID);
354
355   if ( type == SMDSAbs_All )
356     return l.ncells;
357
358   int nb = 0;
359   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
360   for (int i=0; i<l.ncells; i++)
361         {
362            const SMDS_MeshElement* elem = mesh->FindElement(mesh->fromVtkToSmds(l.cells[i]));
363            if (elem->GetType() == type)
364                 nb++;
365         }
366   return nb;
367 }
368
369 ///////////////////////////////////////////////////////////////////////////////
370 /// To be used with STL set
371 ///////////////////////////////////////////////////////////////////////////////
372 bool operator<(const SMDS_MeshNode& e1, const SMDS_MeshNode& e2)
373 {
374         return e1.getVtkId()<e2.getVtkId();
375         /*if(e1.myX<e2.myX) return true;
376         else if(e1.myX==e2.myX)
377         {
378                 if(e1.myY<e2.myY) return true;
379                 else if(e1.myY==e2.myY) return (e1.myZ<e2.myZ);
380                 else return false;
381         }
382         else return false;*/
383 }
384