Salome HOME
Updated copyright comment
[modules/gui.git] / src / VTKViewer / VTKViewer_CellCenters.cxx
1 // Copyright (C) 2007-2024  CEA, EDF, 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 #include "VTKViewer_CellCenters.h"
24
25 #include <vtkCell.h>
26 #include <vtkCellData.h>
27 #include <vtkDataSet.h>
28 #include <vtkInformation.h>
29 #include <vtkInformationVector.h>
30 #include <vtkObjectFactory.h>
31 #include <vtkPointData.h>
32 #include <vtkPoints.h>
33 #include <vtkPolyData.h>
34 #include <vtkCellArray.h>
35
36 vtkStandardNewMacro(VTKViewer_CellCenters)
37
38 /*!
39  * Class       : VTKViewer_CellCenters
40  * Description : Filter computing geometrical centers of given cells
41  *               (differs from native vtk filter by small fix for VTK_CONVEX_POINT_SET cells)
42  */
43
44 /*!
45   Constructor
46 */
47 VTKViewer_CellCenters::VTKViewer_CellCenters()
48 {
49 }
50
51 /*!
52   Redefined main method
53 */
54 int VTKViewer_CellCenters::RequestData(
55   vtkInformation *vtkNotUsed(request),
56   vtkInformationVector **inputVector,
57   vtkInformationVector *outputVector)
58 {
59   // get the info objects
60   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
61   vtkInformation *outInfo = outputVector->GetInformationObject(0);
62
63   // get the input and ouptut
64   vtkDataSet *input = vtkDataSet::SafeDownCast(
65     inInfo->Get(vtkDataObject::DATA_OBJECT()));
66   vtkPolyData *output = vtkPolyData::SafeDownCast(
67     outInfo->Get(vtkDataObject::DATA_OBJECT()));
68
69   vtkIdType cellId, numCells;
70   int subId;
71   vtkCellData *inCD;
72   vtkPointData *outPD;
73   vtkPoints *newPts;
74   vtkCell *cell;
75   double x[3], pcoords[3];
76   double *weights;
77
78   inCD=input->GetCellData();
79   outPD=output->GetPointData();
80
81   if ( (numCells = input->GetNumberOfCells()) < 1 )
82     {
83     vtkWarningMacro(<<"No cells to generate center points for");
84     return 1;
85     }
86
87   newPts = vtkPoints::New();
88   newPts->SetNumberOfPoints(numCells);
89   weights = new double [input->GetMaxCellSize()];
90
91   int abort=0;
92   vtkIdType progressInterval = numCells/10 + 1;
93   int hasEmptyCells = 0;
94   for (cellId=0; cellId < numCells && !abort; cellId++)
95     {
96     if ( ! (cellId % progressInterval) ) 
97       {
98       vtkDebugMacro(<<"Processing #" << cellId);
99       this->UpdateProgress (0.5*cellId/numCells);
100       abort = this->GetAbortExecute();
101       }
102
103     cell = input->GetCell(cellId);
104     if (cell->GetCellType() != VTK_EMPTY_CELL)
105       {
106         // fix for VTK_CONVEX_POINT_SET cells
107         if (cell->GetCellType() == VTK_CONVEX_POINT_SET )
108         {
109           x[0] = x[1] = x[2] = 0;
110           vtkPoints* aPoints = cell->GetPoints();
111           int aNbPoints = aPoints->GetNumberOfPoints();
112           for( int i = 0; i < aNbPoints; i++ )
113           {
114             double aCoord[3];
115             aPoints->GetPoint( i, aCoord );
116             x[0] += aCoord[0];
117             x[1] += aCoord[1];
118             x[2] += aCoord[2];
119           }
120           x[0] /= aNbPoints;
121           x[1] /= aNbPoints;
122           x[2] /= aNbPoints;
123         }
124         else
125         {
126           subId = cell->GetParametricCenter(pcoords);
127           cell->EvaluateLocation(subId, pcoords, x, weights);
128         }
129         newPts->SetPoint(cellId,x);
130       }
131     else
132       {
133       hasEmptyCells = 1;
134       }
135     }
136
137   if ( this->VertexCells )
138     {
139     vtkIdType pts[1];
140     vtkCellData *outCD=output->GetCellData();
141     vtkCellArray *verts = vtkCellArray::New();
142     verts->Allocate(verts->EstimateSize(1,numCells),1);
143
144     for (cellId=0; cellId < numCells && !abort; cellId++)
145       {
146       if ( ! (cellId % progressInterval) ) 
147         {
148         vtkDebugMacro(<<"Processing #" << cellId);
149         this->UpdateProgress (0.5+0.5*cellId/numCells);
150         abort = this->GetAbortExecute();
151         }
152
153       cell = input->GetCell(cellId);
154       if (cell->GetCellType() != VTK_EMPTY_CELL)
155         {
156         pts[0] = cellId;
157         verts->InsertNextCell(1,pts);
158         }
159       }
160
161     output->SetVerts(verts);
162     verts->Delete();
163     if (!hasEmptyCells)
164       {
165       outCD->PassData(inCD); //only if verts are generated
166       }
167     }
168
169   // clean up and update output
170   output->SetPoints(newPts);
171   newPts->Delete();
172
173   if (!hasEmptyCells)
174     {
175     outPD->PassData(inCD); //because number of points = number of cells
176     }
177   if (weights)
178     {
179     delete [] weights;
180     }
181
182   return 1;
183 }