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