Salome HOME
Copyright update 2022
[modules/paravis.git] / src / Plugins / TableReader / plugin / TableTo3DModule / vtkTableTo3D.cxx
1 // Copyright (C) 2010-2022  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "vtkTableTo3D.h"
21
22 #include <vtkContourFilter.h>
23 #include <vtkDoubleArray.h>
24 #include <vtkInformation.h>
25 #include <vtkObjectFactory.h>
26 #include <vtkPointData.h>
27 #include <vtkPoints.h>
28 #include <vtkPolyData.h>
29 #include <vtkSmartPointer.h>
30 #include <vtkStructuredGrid.h>
31 #include <vtkStructuredGridGeometryFilter.h>
32 #include <vtkTable.h>
33 #include <vtkVariantArray.h>
34 #include <vtkWarpScalar.h>
35
36 vtkStandardNewMacro(vtkTableTo3D);
37
38 vtkTableTo3D::vtkTableTo3D()
39 {
40   this->ScaleFactor = 1.0;
41   this->UseOptimusScale = true;
42   this->PresentationType = TABLETO3D_SURFACE;
43   this->NumberOfContours = 32;
44 }
45
46 vtkTableTo3D::~vtkTableTo3D() {}
47
48 int vtkTableTo3D::FillInputPortInformation(int vtkNotUsed(port), vtkInformation* info)
49 {
50   info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
51   return 1;
52 }
53
54 int vtkTableTo3D::RequestData(vtkInformation* vtkNotUsed(request),
55   vtkInformationVector** inputVector, vtkInformationVector* outputVector)
56 {
57   vtkTable* input = vtkTable::GetData(inputVector[0], 0);
58   vtkPolyData* output = vtkPolyData::GetData(outputVector, 0);
59
60   if (input->GetNumberOfRows() == 0 || input->GetNumberOfColumns() < 2)
61   {
62     return 1;
63   }
64
65   vtkIdType xSize = input->GetNumberOfRows();
66   vtkIdType ySize = input->GetNumberOfColumns() - 1;
67   vtkIdType nbPoints = xSize * ySize;
68
69   vtkDataArray* xAxis = vtkDataArray::SafeDownCast(input->GetColumn(0));
70   if (!xAxis)
71   {
72     vtkErrorMacro("The first column is not numeric.");
73     return 1;
74   }
75
76   double xRange = xAxis->GetTuple1(xSize - 1) - xAxis->GetTuple1(0);
77   double yDelta = xRange / ySize;
78
79   vtkSmartPointer<vtkDoubleArray> yAxis = vtkSmartPointer<vtkDoubleArray>::New();
80   yAxis->SetNumberOfValues(ySize);
81   for (vtkIdType i = 0; i < ySize; i++)
82   {
83     yAxis->SetValue(i, i * yDelta);
84   }
85
86   vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
87   points->SetNumberOfPoints(nbPoints);
88
89   vtkSmartPointer<vtkIntArray> pointsIdMapper = vtkSmartPointer<vtkIntArray>::New();
90   pointsIdMapper->SetName("POINTS_ID_MAPPER");
91   pointsIdMapper->SetNumberOfComponents(2);
92   pointsIdMapper->SetNumberOfTuples(nbPoints);
93   int* pointsIdMapperPtr = pointsIdMapper->GetPointer(0);
94
95   for (vtkIdType i = 0, pntId = 0; i < ySize; i++)
96   {
97     for (vtkIdType j = 0; j < xSize; j++, pntId++)
98     {
99       points->SetPoint(pntId, xAxis->GetTuple1(j), yAxis->GetValue(i), 0.0);
100
101       *pointsIdMapperPtr++ = pntId;
102       *pointsIdMapperPtr++ = 0;
103     }
104   }
105
106   vtkSmartPointer<vtkDoubleArray> scalars = vtkSmartPointer<vtkDoubleArray>::New();
107   scalars->SetNumberOfComponents(1);
108   scalars->SetNumberOfTuples(nbPoints);
109   double* scalarsPtr = scalars->GetPointer(0);
110   for (vtkIdType i = 0; i < ySize; i++)
111   {
112     vtkDataArray* col = vtkDataArray::SafeDownCast(input->GetColumn(i + 1));
113
114     if (!col)
115     {
116       vtkErrorMacro("Column " << i << "is not numeric.");
117       return 1;
118     }
119
120     for (vtkIdType j = 0; j < xSize; j++)
121     {
122       double value = col->GetTuple1(j);
123       *scalarsPtr++ = value;
124     }
125   }
126
127   vtkSmartPointer<vtkStructuredGrid> structuredGrid = vtkSmartPointer<vtkStructuredGrid>::New();
128   structuredGrid->SetPoints(points);
129
130   structuredGrid->SetDimensions(xSize, ySize, 1);
131
132   // structuredGrid->GetPointData()->AddArray(pointsIdMapper);
133   if (input->GetInformation()->Has(vtkDataObject::FIELD_NAME()))
134   {
135     scalars->SetName(input->GetInformation()->Get(vtkDataObject::FIELD_NAME()));
136   }
137   else
138   {
139     scalars->SetName("Table");
140   }
141   structuredGrid->GetPointData()->SetScalars(scalars);
142
143   vtkSmartPointer<vtkStructuredGridGeometryFilter> geomFilter =
144     vtkSmartPointer<vtkStructuredGridGeometryFilter>::New();
145   geomFilter->SetInputData(structuredGrid);
146   geomFilter->Update();
147
148   vtkSmartPointer<vtkWarpScalar> warpScalar = vtkSmartPointer<vtkWarpScalar>::New();
149
150   double scaleFactor = this->ScaleFactor;
151   if (this->UseOptimusScale)
152   {
153     double range[2];
154     geomFilter->GetOutput()->GetScalarRange(range);
155     double length = geomFilter->GetOutput()->GetLength();
156     if (range[1] > 0)
157     {
158       scaleFactor = length / range[1] * 0.3;
159     }
160     else
161     {
162       scaleFactor = 0;
163     }
164   }
165
166   if (this->PresentationType == TABLETO3D_SURFACE)
167   {
168     warpScalar->SetInputConnection(geomFilter->GetOutputPort(0));
169     warpScalar->SetScaleFactor(scaleFactor);
170   }
171   else
172   {
173     vtkSmartPointer<vtkContourFilter> contourFilter = vtkSmartPointer<vtkContourFilter>::New();
174     contourFilter->SetInputConnection(geomFilter->GetOutputPort(0));
175     contourFilter->GenerateValues(
176       this->NumberOfContours, geomFilter->GetOutput()->GetScalarRange());
177     warpScalar->SetInputConnection(contourFilter->GetOutputPort(0));
178     warpScalar->SetScaleFactor(scaleFactor);
179   }
180
181   warpScalar->Update();
182   output->ShallowCopy(warpScalar->GetPolyDataOutput());
183
184   return 1;
185 }
186
187 //----------------------------------------------------------------------------
188 void vtkTableTo3D::PrintSelf(ostream& os, vtkIndent indent)
189 {
190   this->Superclass::PrintSelf(os, indent);
191
192   os << indent << "ScaleFactor: " << this->ScaleFactor << endl;
193   os << indent << "UseOptimusScale: " << (this->UseOptimusScale ? "true" : "false") << endl;
194   os << indent << "PresentationType: "
195      << ((this->PresentationType == TABLETO3D_SURFACE) ? "Surface" : "Contour") << endl;
196   os << indent << "NumberOfContours: " << this->NumberOfContours << endl;
197 }