OGS
VtkImageDataToPointCloudFilter.cpp
Go to the documentation of this file.
1 
11 // ** INCLUDES **
13 
14 #include <vtkDoubleArray.h>
15 #include <vtkIdList.h>
16 #include <vtkImageData.h>
17 #include <vtkInformation.h>
18 #include <vtkInformationVector.h>
19 #include <vtkLine.h>
20 #include <vtkObjectFactory.h>
21 #include <vtkPointData.h>
22 #include <vtkPolyData.h>
23 #include <vtkSmartPointer.h>
24 
25 #include <algorithm>
26 #include <vector>
27 
29 
31 
32 void VtkImageDataToPointCloudFilter::PrintSelf(ostream& os, vtkIndent indent)
33 {
34  this->Superclass::PrintSelf(os, indent);
35 }
36 
38  int /*port*/, vtkInformation* info)
39 {
40  info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkImageData");
41  return 1;
42 }
43 
45  vtkInformation* /*request*/,
46  vtkInformationVector** inputVector,
47  vtkInformationVector* outputVector)
48 {
49  vtkDebugMacro(<< "Executing VtkImageDataToPointCloudFilter");
50 
51  vtkInformation* input_info = inputVector[0]->GetInformationObject(0);
52  vtkInformation* output_info = outputVector->GetInformationObject(0);
53  vtkImageData* input = vtkImageData::SafeDownCast(
54  input_info->Get(vtkDataObject::DATA_OBJECT()));
55  vtkPolyData* output = vtkPolyData::SafeDownCast(
56  output_info->Get(vtkDataObject::DATA_OBJECT()));
57 
58  if (input->GetScalarSize() != sizeof(float))
59  {
60  vtkDebugMacro(
61  "Existing data does not have float-type and cannot be processed. "
62  "Aborting filter process...");
63  return 0;
64  }
65  float* pixvals = static_cast<float*>(input->GetScalarPointer());
66  int const n_comp = input->GetNumberOfScalarComponents();
67  if (n_comp < 1)
68  {
69  vtkDebugMacro("Error reading number of components.");
70  }
71  if (n_comp > 2)
72  {
73  vtkDebugMacro(
74  "RGB colours detected. Using only first channel for computation.");
75  }
76 
77  vtkIdType const n_points = input->GetNumberOfPoints();
78  if (n_points == 0)
79  {
80  vtkDebugMacro("No data found!");
81  return -1;
82  }
83 
84  double spacing[3];
85  input->GetSpacing(spacing);
86 
87  if (MinValueRange == -1 || MaxValueRange == -1)
88  {
89  double range[2];
90  vtkPointData* input_data = input->GetPointData();
91  input_data->GetArray(0)->GetRange(range);
92  MinValueRange = range[0];
93  MaxValueRange = range[1];
94  }
95 
96  std::vector<vtkIdType> density;
97  density.reserve(n_points);
98  for (vtkIdType i = 0; i < n_points; ++i)
99  {
100  // Skip transparent pixels
101  if (n_comp == 2 || n_comp == 4)
102  {
103  if (pixvals[(i + 1) * n_comp - 1] < 0.00000001f)
104  {
105  density.push_back(0);
106  continue;
107  }
108  }
109  double const val(pixvals[i * n_comp]);
110  double const calc_gamma = (IsLinear) ? 1 : Gamma;
111  std::size_t const pnts_per_cell =
112  interpolate(MinValueRange, MaxValueRange, val, calc_gamma);
113  density.push_back(static_cast<vtkIdType>(
114  std::floor(pnts_per_cell * GetPointScaleFactor())));
115  }
116 
117  // Allocate the space needed for output objects
118  auto const sum =
119  std::accumulate(density.begin(), density.end(), vtkIdType{0});
120  vtkSmartPointer<vtkPoints> new_points = vtkSmartPointer<vtkPoints>::New();
121  new_points->SetNumberOfPoints(sum);
122  vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
123  cells->Allocate(sum);
124  vtkSmartPointer<vtkDoubleArray> intensity =
125  vtkSmartPointer<vtkDoubleArray>::New();
126  intensity->Allocate(sum);
127  intensity->SetName("Intensity");
128  double const half_cellsize(spacing[0] / 2.0);
129  std::size_t pnt_idx(0);
130  for (std::size_t i = 0; i < static_cast<std::size_t>(n_points); ++i)
131  {
132  if (density[i] == 0)
133  {
134  continue;
135  }
136  double p[3];
137  input->GetPoint(i, p);
138 
139  MathLib::Point3d min_pnt{std::array<double, 3>{
140  {p[0] - half_cellsize, p[1] - half_cellsize, MinHeight}}};
141  MathLib::Point3d max_pnt{std::array<double, 3>{
142  {p[0] + half_cellsize, p[1] + half_cellsize, MaxHeight}}};
143  createPoints(new_points, cells, pnt_idx, density[i], min_pnt, max_pnt);
144  for (vtkIdType j = 0; j < density[i]; ++j)
145  {
146  intensity->InsertValue(pnt_idx + j, pixvals[i * n_comp]);
147  }
148  pnt_idx += density[i];
149  }
150 
151  output->SetPoints(new_points);
152  output->SetVerts(cells);
153  output->GetPointData()->AddArray(intensity);
154  output->GetPointData()->SetActiveAttribute("Intensity",
155  vtkDataSetAttributes::SCALARS);
156  output->Squeeze();
157 
158  vtkDebugMacro(<< "Created " << new_points->GetNumberOfPoints()
159  << " points.");
160  return 1;
161 }
162 
164  vtkSmartPointer<vtkPoints>& points,
165  vtkSmartPointer<vtkCellArray>& cells,
166  std::size_t pnt_idx,
167  std::size_t n_points,
168  MathLib::Point3d const& min_pnt,
169  MathLib::Point3d const& max_pnt) const
170 {
171  for (std::size_t i = 0; i < n_points; ++i)
172  {
173  double p[3];
174  p[0] = getRandomNumber(min_pnt[0], max_pnt[0]);
175  p[1] = getRandomNumber(min_pnt[1], max_pnt[1]);
176  p[2] = getRandomNumber(min_pnt[2], max_pnt[2]);
177  points->SetPoint(pnt_idx + i, p);
178  cells->InsertNextCell(1);
179  cells->InsertCellPoint(pnt_idx + i);
180  }
181 }
182 
184  double const& max) const
185 {
186  return (static_cast<double>(std::rand()) / RAND_MAX) * (max - min) + min;
187 }
188 
190  double high,
191  double p,
192  double gamma) const
193 {
194  p = std::clamp(p, low, high);
195  double const r = (p - low) / (high - low);
196  return static_cast<std::size_t>(
198  std::pow(r, gamma) +
200 }
vtkStandardNewMacro(VtkImageDataToPointCloudFilter)
int FillInputPortInformation(int port, vtkInformation *info) override
Sets input port to vtkImageData.
void PrintSelf(ostream &os, vtkIndent indent) override
Prints information about the filter.
void createPoints(vtkSmartPointer< vtkPoints > &points, vtkSmartPointer< vtkCellArray > &cells, std::size_t pnt_idx, std::size_t n_points, MathLib::Point3d const &min_pnt, MathLib::Point3d const &max_pnt) const
Creates the point objects based on the parameters set by the user.
double getRandomNumber(double const &min, double const &max) const
Returns a random number in [min, max].
std::size_t interpolate(double low, double high, double p, double gamma) const
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Updates the graphical object.
static const double p
static const double r