]> SALOME platform Git repositories - modules/smesh.git/blob - src/SMDS/SMDS_VtkVolume.cxx
Salome HOME
use vtkIdType to be independant of architecture
[modules/smesh.git] / src / SMDS / SMDS_VtkVolume.cxx
1 #include "SMDS_VtkVolume.hxx"
2 #include "SMDS_MeshNode.hxx"
3 #include "SMDS_Mesh.hxx"
4 #include "SMDS_VtkCellIterator.hxx"
5
6 #include "utilities.h"
7
8 #include <vector>
9
10 SMDS_VtkVolume::SMDS_VtkVolume()
11 {
12 }
13
14 SMDS_VtkVolume::SMDS_VtkVolume(std::vector<vtkIdType> nodeIds, SMDS_Mesh* mesh)
15 {
16   init(nodeIds, mesh);
17 }
18 /*!
19  * typed used are vtk types (@see vtkCellType.h)
20  * see GetEntityType() for conversion in SMDS type (@see SMDSAbs_ElementType.hxx)
21  */
22 void SMDS_VtkVolume::init(std::vector<vtkIdType> nodeIds, SMDS_Mesh* mesh)
23 {
24   vtkUnstructuredGrid* grid = mesh->getGrid();
25   myIdInShape = -1;
26   myMeshId = mesh->getMeshId();
27   vtkIdType aType = VTK_TETRA;
28   switch (nodeIds.size())
29   {
30     case 4:
31       aType = VTK_TETRA;
32       break;
33     case 5:
34       aType = VTK_PYRAMID;
35       break;
36     case 6:
37       aType = VTK_WEDGE;
38       break;
39     case 8:
40       aType = VTK_HEXAHEDRON;
41       break;
42     case 10:
43       aType = VTK_QUADRATIC_TETRA;
44       break;
45     case 13:
46       aType = VTK_QUADRATIC_PYRAMID;
47       break;
48     case 15:
49       aType = VTK_QUADRATIC_WEDGE;
50       break;
51     case 20:
52       aType = VTK_QUADRATIC_HEXAHEDRON;
53       break;
54     default:
55       aType = VTK_HEXAHEDRON;
56       break;
57   }
58   myVtkID = grid->InsertNextLinkedCell(aType, nodeIds.size(), &nodeIds[0]);
59   mesh->setMyModified();
60   //MESSAGE("SMDS_VtkVolume::init myVtkID " << myVtkID);
61 }
62
63 //#ifdef VTK_HAVE_POLYHEDRON
64 void SMDS_VtkVolume::initPoly(std::vector<vtkIdType> nodeIds, std::vector<int> nbNodesPerFace, SMDS_Mesh* mesh)
65 {
66   //MESSAGE("SMDS_VtkVolume::initPoly");
67   SMDS_UnstructuredGrid* grid = mesh->getGrid();
68   double center[3];
69   this->gravityCenter(grid, &nodeIds[0], nodeIds.size(), &center[0]);
70   vector<vtkIdType> ptIds;
71   ptIds.clear();
72   vtkIdType nbFaces = nbNodesPerFace.size();
73   int k = 0;
74   for (int i = 0; i < nbFaces; i++)
75     {
76       int nf = nbNodesPerFace[i];
77       ptIds.push_back(nf);
78       double a[3];
79       double b[3];
80       double c[3];
81       grid->GetPoints()->GetPoint(nodeIds[k], a);
82       grid->GetPoints()->GetPoint(nodeIds[k + 1], b);
83       grid->GetPoints()->GetPoint(nodeIds[k + 2], c);
84       bool isFaceForward = this->isForward(a, b, c, center);
85       //MESSAGE("isFaceForward " << i << " " << isFaceForward);
86       vtkIdType *facePts = &nodeIds[k];
87       if (isFaceForward)
88         for (int n = 0; n < nf; n++)
89           ptIds.push_back(facePts[n]);
90       else
91         for (int n = nf - 1; n >= 0; n--)
92           ptIds.push_back(facePts[n]);
93       k += nf;
94     }
95   myVtkID = grid->InsertNextLinkedCell(VTK_POLYHEDRON, nbFaces, &ptIds[0]);
96   mesh->setMyModified();
97 }
98 //#endif
99
100 bool SMDS_VtkVolume::ChangeNodes(const SMDS_MeshNode* nodes[], const int nbNodes)
101 {
102   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
103   vtkIdType npts = 0;
104   vtkIdType* pts = 0;
105   grid->GetCellPoints(myVtkID, npts, pts);
106   if (nbNodes != npts)
107     {
108       MESSAGE("ChangeNodes problem: not the same number of nodes " << npts << " -> " << nbNodes);
109       return false;
110     }
111   for (int i = 0; i < nbNodes; i++)
112     {
113       pts[i] = nodes[i]->getVtkId();
114     }
115   SMDS_Mesh::_meshList[myMeshId]->setMyModified();
116   return true;
117 }
118
119 /*!
120  * Reorder in VTK order a list of nodes given in SMDS order.
121  * To be used before ChangeNodes: lists are given or computed in SMDS order.
122  */
123 bool SMDS_VtkVolume::vtkOrder(const SMDS_MeshNode* nodes[], const int nbNodes)
124 {
125   if (nbNodes != this->NbNodes())
126     {
127       MESSAGE("vtkOrder, wrong number of nodes " << nbNodes << " instead of "<< this->NbNodes());
128       return false;
129     }
130   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
131   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
132   switch (aVtkType)
133   {
134     case VTK_TETRA:
135       this->exchange(nodes, 1, 2);
136       break;
137     case VTK_QUADRATIC_TETRA:
138       this->exchange(nodes, 1, 2);
139       this->exchange(nodes, 4, 6);
140       this->exchange(nodes, 8, 9);
141       break;
142     case VTK_PYRAMID:
143       this->exchange(nodes, 1, 3);
144       break;
145     case VTK_WEDGE:
146       break;
147     case VTK_QUADRATIC_PYRAMID:
148       this->exchange(nodes, 1, 3);
149       this->exchange(nodes, 5, 8);
150       this->exchange(nodes, 6, 7);
151       this->exchange(nodes, 10, 12);
152       break;
153     case VTK_QUADRATIC_WEDGE:
154       break;
155     case VTK_HEXAHEDRON:
156       this->exchange(nodes, 1, 3);
157       this->exchange(nodes, 5, 7);
158       break;
159     case VTK_QUADRATIC_HEXAHEDRON:
160       this->exchange(nodes, 1, 3);
161       this->exchange(nodes, 5, 7);
162       this->exchange(nodes, 8, 11);
163       this->exchange(nodes, 9, 10);
164       this->exchange(nodes, 12, 15);
165       this->exchange(nodes, 13, 14);
166       this->exchange(nodes, 17, 19);
167       break;
168     case VTK_POLYHEDRON:
169     default:
170       break;
171   }
172   return true;
173 }
174
175 SMDS_VtkVolume::~SMDS_VtkVolume()
176 {
177 }
178
179 void SMDS_VtkVolume::Print(ostream & OS) const
180 {
181   OS << "volume <" << GetID() << "> : ";
182 }
183
184 int SMDS_VtkVolume::NbFaces() const
185 {
186   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
187   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
188   int nbFaces = 4;
189   switch (aVtkType)
190   {
191     case VTK_TETRA:
192     case VTK_QUADRATIC_TETRA:
193       nbFaces = 4;
194       break;
195     case VTK_PYRAMID:
196     case VTK_WEDGE:
197     case VTK_QUADRATIC_PYRAMID:
198     case VTK_QUADRATIC_WEDGE:
199       nbFaces = 5;
200       break;
201     case VTK_HEXAHEDRON:
202     case VTK_QUADRATIC_HEXAHEDRON:
203       nbFaces = 6;
204       break;
205     case VTK_POLYHEDRON:
206       {
207         vtkIdType nFaces = 0;
208         vtkIdType* ptIds = 0;
209         grid->GetFaceStream(this->myVtkID, nFaces, ptIds);
210         nbFaces = nFaces;
211         break;
212       }
213     default:
214       MESSAGE("invalid volume type")
215       ;
216       nbFaces = 0;
217       break;
218   }
219   return nbFaces;
220 }
221
222 int SMDS_VtkVolume::NbNodes() const
223 {
224   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
225   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
226   int nbPoints = 0;
227   if (aVtkType != VTK_POLYHEDRON)
228     {
229       nbPoints = grid->GetCell(myVtkID)->GetNumberOfPoints();
230     }
231   else
232     {
233       vtkIdType nFaces = 0;
234       vtkIdType* ptIds = 0;
235       grid->GetFaceStream(this->myVtkID, nFaces, ptIds);
236       int id = 0;
237       for (int i = 0; i < nFaces; i++)
238         {
239           int nodesInFace = ptIds[id];
240           nbPoints += nodesInFace;
241           id += (nodesInFace + 1);
242         }
243     }
244   return nbPoints;
245 }
246
247 int SMDS_VtkVolume::NbEdges() const
248 {
249   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
250   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
251   int nbEdges = 6;
252   switch (aVtkType)
253   {
254     case VTK_TETRA:
255     case VTK_QUADRATIC_TETRA:
256       nbEdges = 6;
257       break;
258     case VTK_PYRAMID:
259     case VTK_QUADRATIC_PYRAMID:
260       nbEdges = 8;
261       break;
262     case VTK_WEDGE:
263     case VTK_QUADRATIC_WEDGE:
264       nbEdges = 9;
265       break;
266     case VTK_HEXAHEDRON:
267     case VTK_QUADRATIC_HEXAHEDRON:
268       nbEdges = 12;
269       break;
270     case VTK_POLYHEDRON:
271       {
272         vtkIdType nFaces = 0;
273         vtkIdType* ptIds = 0;
274         grid->GetFaceStream(this->myVtkID, nFaces, ptIds);
275         nbEdges = 0;
276         int id = 0;
277         for (int i = 0; i < nFaces; i++)
278           {
279             int edgesInFace = ptIds[id];
280             id += (edgesInFace + 1);
281             nbEdges += edgesInFace;
282           }
283         nbEdges = nbEdges / 2;
284         break;
285       }
286     default:
287       MESSAGE("invalid volume type")
288       ;
289       nbEdges = 0;
290       break;
291   }
292   return nbEdges;
293 }
294
295 /*! polyhedron only,
296  *  1 <= face_ind <= NbFaces()
297  */
298 int SMDS_VtkVolume::NbFaceNodes(const int face_ind) const
299 {
300   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
301   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
302   int nbNodes = 0;
303   if (aVtkType == VTK_POLYHEDRON)
304     {
305       vtkIdType nFaces = 0;
306       vtkIdType* ptIds = 0;
307       grid->GetFaceStream(this->myVtkID, nFaces, ptIds);
308       int id = 0;
309       for (int i = 0; i < nFaces; i++)
310         {
311           int nodesInFace = ptIds[id];
312           id += (nodesInFace + 1);
313           if (i == face_ind - 1)
314             {
315               nbNodes = nodesInFace;
316               break;
317             }
318         }
319     }
320   return nbNodes;
321 }
322
323 /*! polyhedron only,
324  *  1 <= face_ind <= NbFaces()
325  *  1 <= node_ind <= NbFaceNodes()
326  */
327 const SMDS_MeshNode* SMDS_VtkVolume::GetFaceNode(const int face_ind, const int node_ind) const
328 {
329   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
330   vtkUnstructuredGrid* grid = mesh->getGrid();
331   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
332   const SMDS_MeshNode* node = 0;
333   if (aVtkType == VTK_POLYHEDRON)
334     {
335       vtkIdType nFaces = 0;
336       vtkIdType* ptIds = 0;
337       grid->GetFaceStream(this->myVtkID, nFaces, ptIds);
338       int id = 0;
339       for (int i = 0; i < nFaces; i++)
340         {
341           int nodesInFace = ptIds[id]; // nodeIds in ptIds[id+1 .. id+nodesInFace]
342           if (i == face_ind - 1) // first face is number 1
343             {
344               if ((node_ind > 0) && (node_ind <= nodesInFace))
345                 node = mesh->FindNodeVtk(ptIds[id + node_ind]); // ptIds[id+1] : first node
346               break;
347             }
348           id += (nodesInFace + 1);
349         }
350     }
351   return node;
352 }
353
354 /*! polyhedron only,
355  *  return number of nodes for each face
356  */
357 const std::vector<int> SMDS_VtkVolume::GetQuantities() const
358 {
359   vector<int> quantities;
360   quantities.clear();
361   SMDS_Mesh *mesh = SMDS_Mesh::_meshList[myMeshId];
362   vtkUnstructuredGrid* grid = mesh->getGrid();
363   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
364   if (aVtkType == VTK_POLYHEDRON)
365     {
366       vtkIdType nFaces = 0;
367       vtkIdType* ptIds = 0;
368       grid->GetFaceStream(this->myVtkID, nFaces, ptIds);
369       int id = 0;
370       for (int i = 0; i < nFaces; i++)
371         {
372           int nodesInFace = ptIds[id]; // nodeIds in ptIds[id+1 .. id+nodesInFace]
373           quantities.push_back(nodesInFace);
374           id += (nodesInFace + 1);
375         }
376     }
377   return quantities;
378 }
379
380 SMDS_ElemIteratorPtr SMDS_VtkVolume::elementsIterator(SMDSAbs_ElementType type) const
381 {
382   switch (type)
383   {
384     case SMDSAbs_Node:
385       {
386         SMDSAbs_EntityType aType = this->GetEntityType();
387         if (aType == SMDSEntity_Polyhedra)
388           return SMDS_ElemIteratorPtr(new SMDS_VtkCellIteratorPolyH(SMDS_Mesh::_meshList[myMeshId], myVtkID, aType));
389         else
390           return SMDS_ElemIteratorPtr(new SMDS_VtkCellIterator(SMDS_Mesh::_meshList[myMeshId], myVtkID, aType));
391       }
392     default:
393       MESSAGE("ERROR : Iterator not implemented");
394       return SMDS_ElemIteratorPtr((SMDS_ElemIterator*) NULL);
395   }
396 }
397
398 SMDS_ElemIteratorPtr SMDS_VtkVolume::nodesIteratorToUNV() const
399 {
400   return SMDS_ElemIteratorPtr(new SMDS_VtkCellIteratorToUNV(SMDS_Mesh::_meshList[myMeshId], myVtkID, GetEntityType()));
401 }
402
403 SMDS_ElemIteratorPtr SMDS_VtkVolume::interlacedNodesElemIterator() const
404 {
405   return SMDS_ElemIteratorPtr(new SMDS_VtkCellIteratorToUNV(SMDS_Mesh::_meshList[myMeshId], myVtkID, GetEntityType()));
406 }
407
408 SMDSAbs_ElementType SMDS_VtkVolume::GetType() const
409 {
410   return SMDSAbs_Volume;
411 }
412
413 /*!
414  * \brief Return node by its index
415  * \param ind - node index
416  * \retval const SMDS_MeshNode* - the node
417  */
418 const SMDS_MeshNode* SMDS_VtkVolume::GetNode(const int ind) const
419 {
420   // TODO optimize if possible (vtkCellIterator)
421   return SMDS_MeshElement::GetNode(ind);
422 }
423
424 bool SMDS_VtkVolume::IsQuadratic() const
425 {
426   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
427   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
428   // TODO quadratic polyhedrons ?
429   switch (aVtkType)
430   {
431     case VTK_QUADRATIC_TETRA:
432     case VTK_QUADRATIC_PYRAMID:
433     case VTK_QUADRATIC_WEDGE:
434     case VTK_QUADRATIC_HEXAHEDRON:
435       return true;
436       break;
437     default:
438       return false;
439   }
440 }
441
442 bool SMDS_VtkVolume::IsPoly() const
443 {
444   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
445   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
446   return (aVtkType == VTK_POLYHEDRON);
447 }
448
449 bool SMDS_VtkVolume::IsMediumNode(const SMDS_MeshNode* node) const
450 {
451   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
452   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
453   int rankFirstMedium = 0;
454   switch (aVtkType)
455   {
456     case VTK_QUADRATIC_TETRA:
457       rankFirstMedium = 4; // medium nodes are of rank 4 to 9
458       break;
459     case VTK_QUADRATIC_PYRAMID:
460       rankFirstMedium = 5; // medium nodes are of rank 5 to 12
461       break;
462     case VTK_QUADRATIC_WEDGE:
463       rankFirstMedium = 6; // medium nodes are of rank 6 to 14
464       break;
465     case VTK_QUADRATIC_HEXAHEDRON:
466       rankFirstMedium = 8; // medium nodes are of rank 8 to 19
467       break;
468     default:
469       return false;
470   }
471   vtkIdType npts = 0;
472   vtkIdType* pts = 0;
473   grid->GetCellPoints(myVtkID, npts, pts);
474   vtkIdType nodeId = node->getVtkId();
475   for (int rank = 0; rank < npts; rank++)
476     {
477       if (pts[rank] == nodeId)
478         {
479           if (rank < rankFirstMedium)
480             return false;
481           else
482             return true;
483         }
484     }
485   //throw SALOME_Exception(LOCALIZED("node does not belong to this element"));
486   MESSAGE("======================================================");
487   MESSAGE("= IsMediumNode: node does not belong to this element =");
488   MESSAGE("======================================================");
489   return false;
490 }
491
492 SMDSAbs_EntityType SMDS_VtkVolume::GetEntityType() const
493 {
494   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
495   vtkIdType aVtkType = grid->GetCellType(this->myVtkID);
496
497   SMDSAbs_EntityType aType = SMDSEntity_Tetra;
498   switch (aVtkType)
499   {
500     case VTK_TETRA:
501       aType = SMDSEntity_Tetra;
502       break;
503     case VTK_PYRAMID:
504       aType = SMDSEntity_Pyramid;
505       break;
506     case VTK_WEDGE:
507       aType = SMDSEntity_Penta;
508       break;
509     case VTK_HEXAHEDRON:
510       aType = SMDSEntity_Hexa;
511       break;
512     case VTK_QUADRATIC_TETRA:
513       aType = SMDSEntity_Quad_Tetra;
514       break;
515     case VTK_QUADRATIC_PYRAMID:
516       aType = SMDSEntity_Quad_Pyramid;
517       break;
518     case VTK_QUADRATIC_WEDGE:
519       aType = SMDSEntity_Quad_Penta;
520       break;
521     case VTK_QUADRATIC_HEXAHEDRON:
522       aType = SMDSEntity_Quad_Hexa;
523       break;
524 //#ifdef VTK_HAVE_POLYHEDRON
525     case VTK_POLYHEDRON:
526       aType = SMDSEntity_Polyhedra;
527       break;
528 //#endif
529     default:
530       aType = SMDSEntity_Polyhedra;
531       break;
532   }
533   return aType;
534 }
535
536 vtkIdType SMDS_VtkVolume::GetVtkType() const
537 {
538   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
539   vtkIdType aType = grid->GetCellType(myVtkID);
540   return aType;
541 }
542
543 void SMDS_VtkVolume::gravityCenter(SMDS_UnstructuredGrid* grid, vtkIdType *nodeIds, int nbNodes, double* result)
544 {
545   for (int j = 0; j < 3; j++)
546     result[j] = 0;
547   if (nbNodes <= 0)
548     return;
549   for (int i = 0; i < nbNodes; i++)
550     {
551       double *coords = grid->GetPoint(nodeIds[i]);
552       for (int j = 0; j < 3; j++)
553         result[j] += coords[j];
554     }
555   for (int j = 0; j < 3; j++)
556     result[j] = result[j] / nbNodes;
557   //MESSAGE("center " << result[0] << " " << result[1] << " "  << result[2]);
558   return;
559 }
560
561 bool SMDS_VtkVolume::isForward(double* a, double* b, double* c, double* d)
562 {
563   double u[3], v[3], w[3];
564   for (int j = 0; j < 3; j++)
565     {
566       //MESSAGE("a,b,c,d " << a[j] << " " << b[j] << " " << c[j] << " " << d[j]);
567       u[j] = b[j] - a[j];
568       v[j] = c[j] - a[j];
569       w[j] = d[j] - a[j];
570       //MESSAGE("u,v,w " << u[j] << " " << v[j] << " " << w[j]);
571     }
572   double prodmixte = (u[1]*v[2] - u[2]*v[1]) * w[0]
573                    + (u[2]*v[0] - u[0]*v[2]) * w[1]
574                    + (u[0]*v[1] - u[1]*v[0]) * w[2];
575   return (prodmixte < 0);
576 }
577
578 /*! For polyhedron only
579  *  @return actual number of nodes (not the sum of nodes of all faces)
580  */
581 int SMDS_VtkVolume::NbUniqueNodes() const
582 {
583   vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
584   return grid->GetCell(myVtkID)->GetNumberOfPoints();
585 }
586
587 /*! For polyhedron use only
588  *  @return iterator on actual nodes (not through the faces)
589  */
590 SMDS_ElemIteratorPtr SMDS_VtkVolume::uniqueNodesIterator() const
591 {
592   MESSAGE("uniqueNodesIterator");
593   return SMDS_ElemIteratorPtr(new SMDS_VtkCellIterator(SMDS_Mesh::_meshList[myMeshId], myVtkID, GetEntityType()));
594 }