TTK
Loading...
Searching...
No Matches
ttkContinuousScatterPlot.cpp
Go to the documentation of this file.
2#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkCharArray.h>
6#include <vtkDataSet.h>
7#include <vtkDoubleArray.h>
8#include <vtkInformation.h>
9#include <vtkObjectFactory.h>
10#include <vtkPointData.h>
11#include <vtkUnstructuredGrid.h>
12
13using namespace std;
14using namespace ttk;
15
17
19 SetNumberOfInputPorts(1);
20 SetNumberOfOutputPorts(1);
21}
22
24
26 vtkInformation *info) {
27
28 if(port == 0) {
29 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkDataSet");
30 }
31 return 1;
32}
33
35 vtkInformation *info) {
36
37 if(port == 0) {
38 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
39 }
40 return 1;
41}
42
43template <typename dataType1, class triangulationType>
44int ttkContinuousScatterPlot::dispatch(const dataType1 *scalars1,
45 vtkDataArray *inputScalars2,
46 const triangulationType *triangulation) {
47 int status = 0;
48 switch(inputScalars2->GetDataType()) {
49 vtkTemplateMacro(
50 (status = this->execute<dataType1, VTK_TT, triangulationType>(
51 scalars1, (VTK_TT *)ttkUtils::GetVoidPointer(inputScalars2),
52 triangulation)));
53 };
54 return status;
55}
56
58 vtkInformationVector **inputVector,
59 vtkInformationVector *outputVector) {
60 vtkDataSet *input = vtkDataSet::GetData(inputVector[0], 0);
61 vtkDataSet *output = vtkDataSet::GetData(outputVector, 0);
62 if(!input || !output)
63 return 0;
64
66 if(!triangulation)
67 return 0;
68
69 auto inputScalars1 = this->GetInputArrayToProcess(0, input);
70 auto inputScalars2 = this->GetInputArrayToProcess(1, input);
71
72#ifndef TTK_ENABLE_KAMIKAZE
73 // wrong scalar fields
74 if(!inputScalars1 || !inputScalars2) {
75 this->printErr("wrong scalar fields.");
76 return -2;
77 }
78#endif
79
80 SimplexId const numberOfPixels
81 = ScatterplotResolution[0] * ScatterplotResolution[1];
82#ifndef TTK_ENABLE_KAMIKAZE
83 // no pixels
84 if(!numberOfPixels) {
85 this->printErr("no pixels.");
86 return -3;
87 }
88#endif
89
90 std::vector<std::vector<double>> density(ScatterplotResolution[0]);
91 std::vector<std::vector<char>> validPointMask(ScatterplotResolution[0]);
92 for(SimplexId k = 0; k < ScatterplotResolution[0]; ++k) {
93 density[k].resize(ScatterplotResolution[1], 0.0);
94 validPointMask[k].resize(ScatterplotResolution[1], 0);
95 }
96
97 SimplexId const numberOfPoints = input->GetNumberOfPoints();
98#ifndef TTK_ENABLE_KAMIKAZE
99 // no points
100 if(numberOfPoints < 1) {
101 this->printErr("no points.");
102 return -4;
103 }
104#endif
105
106 std::array<double, 2> scalarMin{0, 0};
107 std::array<double, 2> scalarMax{0, 0};
108 for(SimplexId k = 0; k < numberOfPoints; ++k) {
109 double const d1 = inputScalars1->GetTuple1(k);
110 double const d2 = inputScalars2->GetTuple1(k);
111
112 if(!k or scalarMin[0] > d1)
113 scalarMin[0] = d1;
114 if(!k or scalarMin[1] > d2)
115 scalarMin[1] = d2;
116 if(!k or scalarMax[0] < d1)
117 scalarMax[0] = d1;
118 if(!k or scalarMax[1] < d2)
119 scalarMax[1] = d2;
120 }
121#ifndef TTK_ENABLE_KAMIKAZE
122 // scalar fields stats problem
123 if(scalarMin[0] == scalarMax[0] or scalarMin[1] == scalarMax[1]) {
124 this->printErr("scalar fields stats problem.");
125 return -5;
126 }
127#endif
128
129 // calling the executing package
130 this->setVertexNumber(numberOfPoints);
131 this->setDummyValue(WithDummyValue, DummyValue);
132 this->setResolutions(ScatterplotResolution[0], ScatterplotResolution[1]);
133 this->setScalarMin(scalarMin);
134 this->setScalarMax(scalarMax);
135 this->setOutputDensity(&density);
136 this->setOutputMask(&validPointMask);
137
138 int status = 0;
139 ttkVtkTemplateMacro(inputScalars1->GetDataType(), triangulation->getType(),
140 (status = this->dispatch<VTK_TT, TTK_TT>(
141 (VTK_TT *)ttkUtils::GetVoidPointer(inputScalars1),
142 inputScalars2, (TTK_TT *)triangulation->getData())));
143
144 // something wrong in baseCode
145 if(status != 0) {
146 std::stringstream msg;
147 msg << "ContinuousScatterPlot.execute() error " << status;
148 this->printErr(msg.str());
149 return -6;
150 }
151
152 vtkNew<vtkCharArray> maskScalars;
153 maskScalars->SetNumberOfComponents(1);
154 maskScalars->SetNumberOfTuples(numberOfPixels);
155 maskScalars->SetName("ValidPointMask");
156
157 vtkNew<vtkDoubleArray> densityScalars;
158 densityScalars->SetNumberOfComponents(1);
159 densityScalars->SetNumberOfTuples(numberOfPixels);
160 densityScalars->SetName("Density");
161
162 vtkNew<vtkDoubleArray> scalars1;
163 scalars1->SetNumberOfComponents(1);
164 scalars1->SetNumberOfTuples(numberOfPixels);
165 scalars1->SetName(inputScalars1->GetName());
166
167 vtkNew<vtkDoubleArray> scalars2;
168 scalars2->SetNumberOfComponents(1);
169 scalars2->SetNumberOfTuples(numberOfPixels);
170 scalars2->SetName(inputScalars2->GetName());
171
172 double delta[2];
173 delta[0] = (scalarMax[0] - scalarMin[0]) / (ScatterplotResolution[0] - 1);
174 delta[1] = (scalarMax[1] - scalarMin[1]) / (ScatterplotResolution[1] - 1);
175
176 double imageMin[2]{0.0, 0.0};
177 double imageMax[2]{1.0, 1.0};
178 if(ProjectImageSupport) {
179 imageMin[0] = scalarMin[0];
180 imageMin[1] = scalarMin[1];
181 imageMax[0] = scalarMax[0];
182 imageMax[1] = scalarMax[1];
183 }
184 double imageDelta[2];
185 imageDelta[0] = (imageMax[0] - imageMin[0]) / (ScatterplotResolution[0] - 1);
186 imageDelta[1] = (imageMax[1] - imageMin[1]) / (ScatterplotResolution[1] - 1);
187
188 vtkNew<vtkUnstructuredGrid> vtu;
189 vtkNew<vtkPoints> pts;
190 pts->SetNumberOfPoints(numberOfPixels);
191
192 SimplexId id{};
193 vtkIdType ids[3];
194 for(SimplexId i = 0; i < ScatterplotResolution[0]; i++) {
195 for(SimplexId j = 0; j < ScatterplotResolution[1]; j++) {
196 // positions:
197 double const x = imageMin[0] + i * imageDelta[0];
198 double const y = imageMin[1] + j * imageDelta[1];
199 pts->SetPoint(id, x, y, 0);
200
201 // scalars:
202 // valid mask
203 maskScalars->SetTuple1(id, validPointMask[i][j]);
204 // density
205 densityScalars->SetTuple1(id, density[i][j]);
206 // original scalar fields
207 double const d1 = scalarMin[0] + i * delta[0];
208 double const d2 = scalarMin[1] + j * delta[1];
209 scalars1->SetTuple1(id, d1);
210 scalars2->SetTuple1(id, d2);
211 if(i < ScatterplotResolution[0] - 1
212 and j < ScatterplotResolution[1] - 1) {
213 ids[0] = id;
214 ids[1] = id + 1;
215 ids[2] = id + ScatterplotResolution[1];
216 vtu->InsertNextCell(VTK_TRIANGLE, 3, ids);
217
218 ids[0] = id + 1;
219 ids[1] = id + ScatterplotResolution[1];
220 ids[2] = id + ScatterplotResolution[1] + 1;
221 vtu->InsertNextCell(VTK_TRIANGLE, 3, ids);
222 }
223
224 ++id;
225 }
226 }
227 vtu->SetPoints(pts);
228 vtu->GetPointData()->AddArray(maskScalars);
229 vtu->GetPointData()->AddArray(densityScalars);
230 vtu->GetPointData()->AddArray(scalars1);
231 vtu->GetPointData()->AddArray(scalars2);
232 output->ShallowCopy(vtu);
233
234 return 1;
235}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter that computes the continuous scatterplot of bivariate volumetric data.
~ttkContinuousScatterPlot() override
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int setDummyValue(bool withDummyValue, double dummyValue)
int setOutputDensity(std::vector< std::vector< double > > *density)
int setOutputMask(std::vector< std::vector< char > > *mask)
int setResolutions(const SimplexId &resolutionX, const SimplexId &resolutionY)
int setVertexNumber(const SimplexId &vertexNumber)
int setScalarMax(const std::array< double, 2 > &scalarMax)
int setScalarMin(const std::array< double, 2 > &scalarMin)
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkContinuousScatterPlot)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69