Salome HOME
Merge from V6_main_20120808 08Aug12
[modules/paravis.git] / src / Plugins / TableReader / TableTo3DFilter / vtkTableTo3D.cxx
1 // Copyright (C) 2010-2012  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.
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 "vtkSmartPointer.h"
23 #include "vtkDoubleArray.h"
24 #include "vtkVariantArray.h"
25 #include "vtkObjectFactory.h"
26 #include "vtkPointData.h"
27 #include "vtkPoints.h"
28 #include "vtkPolyData.h"
29 #include "vtkTable.h"
30 #include "vtkInformation.h"
31 #include "vtkStructuredGrid.h"
32 #include "vtkStructuredGridGeometryFilter.h"
33 #include "vtkWarpScalar.h"
34 #include "vtkContourFilter.h"
35
36 vtkStandardNewMacro(vtkTableTo3D);
37 vtkCxxRevisionMacro(vtkTableTo3D, "$Revision$");
38
39
40 vtkTableTo3D::vtkTableTo3D()
41 {
42   this->ScaleFactor = 1.0;
43   this->UseOptimusScale = true;
44   this->PresentationType = TABLETO3D_SURFACE;
45   this->NumberOfContours = 32;
46 }
47
48 vtkTableTo3D::~vtkTableTo3D()
49 {
50 }
51
52 int vtkTableTo3D::FillInputPortInformation(
53   int vtkNotUsed(port), vtkInformation* info)
54 {
55   info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
56   return 1;
57 }
58
59 int vtkTableTo3D::RequestData(vtkInformation* vtkNotUsed(request),
60   vtkInformationVector** inputVector, vtkInformationVector* outputVector)
61 {
62   vtkTable* input = vtkTable::GetData(inputVector[0], 0);
63   vtkPolyData* output = vtkPolyData::GetData(outputVector, 0);
64
65   if (input->GetNumberOfRows() == 0 ||input->GetNumberOfColumns() < 2)
66     {
67       return 1;
68     }
69   
70   vtkIdType xSize = input->GetNumberOfRows();
71   vtkIdType ySize = input->GetNumberOfColumns() - 1; 
72   vtkIdType nbPoints = xSize * ySize;
73
74   vtkDataArray* xAxis = vtkDataArray::SafeDownCast(input->GetColumn(0));
75   if (!xAxis)
76     {
77       vtkErrorMacro("The first column is not numeric.");
78       return 1;
79     }
80     
81   double xRange = xAxis->GetTuple1(xSize - 1) - xAxis->GetTuple1(0);
82   double yDelta = xRange / ySize;
83   
84   vtkSmartPointer<vtkDoubleArray> yAxis = 
85     vtkSmartPointer<vtkDoubleArray>::New();
86   yAxis->SetNumberOfValues(ySize);
87   for (vtkIdType i = 0; i < ySize; i++ )
88     {
89       yAxis->SetValue(i, i*yDelta);
90     }
91
92   vtkSmartPointer<vtkPoints> points = 
93     vtkSmartPointer<vtkPoints>::New();
94   points->SetNumberOfPoints(nbPoints);
95
96   vtkSmartPointer<vtkIntArray> pointsIdMapper = 
97     vtkSmartPointer<vtkIntArray>::New();
98   pointsIdMapper->SetName("POINTS_ID_MAPPER");
99   pointsIdMapper->SetNumberOfComponents(2);
100   pointsIdMapper->SetNumberOfTuples(nbPoints);
101   int *pointsIdMapperPtr = pointsIdMapper->GetPointer(0);
102
103   for (vtkIdType i = 0, pntId = 0; i < ySize; i++) 
104     {
105       for (vtkIdType j = 0; j < xSize; j++, pntId++) 
106         {
107           points->SetPoint(pntId, xAxis->GetTuple1(j), 
108                                   yAxis->GetValue(i), 
109                                   0.0);
110
111           *pointsIdMapperPtr++ = pntId;
112           *pointsIdMapperPtr++ = 0;
113         }
114     }
115
116   vtkSmartPointer<vtkDoubleArray> scalars = 
117     vtkSmartPointer<vtkDoubleArray>::New();
118   scalars->SetNumberOfComponents(1);
119   scalars->SetNumberOfTuples(nbPoints);
120   double *scalarsPtr = scalars->GetPointer(0);
121   for (vtkIdType i = 0; i < ySize; i++) 
122     {
123       vtkDataArray* col = 
124         vtkDataArray::SafeDownCast(input->GetColumn(i + 1));
125       
126       if (!col)
127         {
128           vtkErrorMacro("Column "<< i <<"is not numeric.");
129           return 1;
130         }
131       
132       for ( vtkIdType j = 0; j < xSize; j++ ) 
133         {
134           double value = col->GetTuple1(j);
135           *scalarsPtr++ = value;
136         }
137     }
138
139   vtkSmartPointer<vtkStructuredGrid> structuredGrid = 
140     vtkSmartPointer<vtkStructuredGrid>::New();
141   structuredGrid->SetPoints(points);
142
143   structuredGrid->SetDimensions(xSize, ySize, 1);
144
145   // structuredGrid->GetPointData()->AddArray(pointsIdMapper);
146   if (input->GetInformation()->Has(vtkDataObject::FIELD_NAME()))
147     {
148       scalars->SetName(input->GetInformation()->Get(vtkDataObject::FIELD_NAME()));
149     }
150   else
151     {
152       scalars->SetName("Table");
153     }
154   structuredGrid->GetPointData()->SetScalars(scalars);
155
156   vtkSmartPointer<vtkStructuredGridGeometryFilter> geomFilter = 
157     vtkSmartPointer<vtkStructuredGridGeometryFilter>::New();
158   geomFilter->SetInput(structuredGrid);
159   geomFilter->Update();
160   
161   vtkSmartPointer<vtkWarpScalar> warpScalar = 
162     vtkSmartPointer<vtkWarpScalar>::New();
163   
164   double scaleFactor = this->ScaleFactor;
165   if (this->UseOptimusScale)
166     {
167       double range[2];
168       geomFilter->GetOutput()->GetScalarRange(range);
169       double length = geomFilter->GetOutput()->GetLength();
170       if (range[1] > 0)
171         {
172           scaleFactor = length / range[1] * 0.3;
173         }
174       else
175         {
176           scaleFactor = 0;
177         }
178     }
179
180   if (this->PresentationType == TABLETO3D_SURFACE)
181     {
182       warpScalar->SetInput(geomFilter->GetOutput());
183       warpScalar->SetScaleFactor(scaleFactor);
184     }
185   else
186     {
187       vtkSmartPointer<vtkContourFilter> contourFilter = 
188         vtkSmartPointer<vtkContourFilter>::New();
189       contourFilter->SetInput(geomFilter->GetOutput());
190       contourFilter->GenerateValues(this->NumberOfContours, 
191                                     geomFilter->GetOutput()->GetScalarRange());
192       warpScalar->SetInput(contourFilter->GetOutput());
193       warpScalar->SetScaleFactor(scaleFactor);
194     }
195
196   warpScalar->Update();
197   output->ShallowCopy(warpScalar->GetPolyDataOutput());
198   
199   return 1;
200 }
201
202 //----------------------------------------------------------------------------
203 void vtkTableTo3D::PrintSelf(ostream& os, vtkIndent indent)
204 {
205   this->Superclass::PrintSelf(os, indent);
206
207   os << indent << "ScaleFactor: " << this->ScaleFactor << endl;
208   os << indent << "UseOptimusScale: "
209      << (this->UseOptimusScale? "true" : "false") << endl;
210   os << indent << "PresentationType: " 
211      << ((this->PresentationType == TABLETO3D_SURFACE)? "Surface" : "Contour")
212      << endl;
213   os << indent << "NumberOfContours: " << this->NumberOfContours << endl;
214 }