OGS
VtkImageDataToPointCloudFilter Class Reference

Detailed Description

A VTK Filter that will create a point cloud with local densities based on pixel values.

Definition at line 16 of file VtkImageDataToPointCloudFilter.h.

#include <VtkImageDataToPointCloudFilter.h>

Inheritance diagram for VtkImageDataToPointCloudFilter:
[legend]
Collaboration diagram for VtkImageDataToPointCloudFilter:
[legend]

Public Member Functions

 vtkTypeMacro (VtkImageDataToPointCloudFilter, vtkPolyDataAlgorithm)
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints information about the filter.
void useLinearInterpolation ()
void useLogarithmicInterpolation (double gamma)
 vtkGetMacro (IsLinear, bool)
 vtkSetMacro (IsLinear, bool)
 vtkGetMacro (MinNumberOfPointsPerCell, vtkIdType)
 vtkSetMacro (MinNumberOfPointsPerCell, vtkIdType)
 vtkGetMacro (MaxNumberOfPointsPerCell, vtkIdType)
 vtkSetMacro (MaxNumberOfPointsPerCell, vtkIdType)
 vtkGetMacro (Gamma, double)
 vtkSetMacro (Gamma, double)
 vtkGetMacro (PointScaleFactor, double)
 vtkSetMacro (PointScaleFactor, double)
 vtkGetMacro (MinValueRange, double)
 vtkSetMacro (MinValueRange, double)
 vtkGetMacro (MaxValueRange, double)
 vtkSetMacro (MaxValueRange, double)
 vtkGetMacro (MinHeight, double)
 vtkSetMacro (MinHeight, double)
 vtkGetMacro (MaxHeight, double)
 vtkSetMacro (MaxHeight, double)

Static Public Member Functions

static VtkImageDataToPointCloudFilterNew ()
 Create a new objects (required because of VTKs reference counting)

Protected Member Functions

 VtkImageDataToPointCloudFilter ()
 ~VtkImageDataToPointCloudFilter () override=default
int FillInputPortInformation (int port, vtkInformation *info) override
 Sets input port to vtkImageData.
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 Updates the graphical object.

Protected Attributes

double Gamma {1.0}
double PointScaleFactor {1.0}
double MinValueRange {-1}
double MaxValueRange {-1}
double MinHeight {0}
double MaxHeight {1000}
vtkIdType MinNumberOfPointsPerCell {1}
vtkIdType MaxNumberOfPointsPerCell {20}
bool IsLinear {true}

Private Member Functions

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

Constructor & Destructor Documentation

◆ VtkImageDataToPointCloudFilter()

VtkImageDataToPointCloudFilter::VtkImageDataToPointCloudFilter ( )
protecteddefault

Referenced by New(), and vtkTypeMacro().

◆ ~VtkImageDataToPointCloudFilter()

VtkImageDataToPointCloudFilter::~VtkImageDataToPointCloudFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ createPoints()

void VtkImageDataToPointCloudFilter::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
private

Creates the point objects based on the parameters set by the user.

Definition at line 156 of file VtkImageDataToPointCloudFilter.cpp.

163{
164 for (std::size_t i = 0; i < n_points; ++i)
165 {
166 double p[3];
167 p[0] = getRandomNumber(min_pnt[0], max_pnt[0]);
168 p[1] = getRandomNumber(min_pnt[1], max_pnt[1]);
169 p[2] = getRandomNumber(min_pnt[2], max_pnt[2]);
170 points->SetPoint(pnt_idx + i, p);
171 cells->InsertNextCell(1);
172 cells->InsertCellPoint(pnt_idx + i);
173 }
174}
double getRandomNumber(double const &min, double const &max) const
Returns a random number in [min, max].

References getRandomNumber().

Referenced by RequestData().

◆ FillInputPortInformation()

int VtkImageDataToPointCloudFilter::FillInputPortInformation ( int port,
vtkInformation * info )
overrideprotected

Sets input port to vtkImageData.

Definition at line 30 of file VtkImageDataToPointCloudFilter.cpp.

32{
33 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkImageData");
34 return 1;
35}

◆ getRandomNumber()

double VtkImageDataToPointCloudFilter::getRandomNumber ( double const & min,
double const & max ) const
private

Returns a random number in [min, max].

Definition at line 176 of file VtkImageDataToPointCloudFilter.cpp.

178{
179 return (static_cast<double>(std::rand()) / RAND_MAX) * (max - min) + min;
180}

Referenced by createPoints().

◆ interpolate()

std::size_t VtkImageDataToPointCloudFilter::interpolate ( double low,
double high,
double p,
double gamma ) const
private

Interpolates the required number of points. Using gamma = 1 this is a linear interpolation, for values smaller/greater than 1 the curve becomes logarithmic/exponential, resulting in a smaller/larger difference of point densities given the same input values.

Definition at line 182 of file VtkImageDataToPointCloudFilter.cpp.

186{
187 p = std::clamp(p, low, high);
188 double const r = (p - low) / (high - low);
189 return static_cast<std::size_t>(
191 std::pow(r, gamma) +
193}

References MaxNumberOfPointsPerCell, and MinNumberOfPointsPerCell.

Referenced by RequestData().

◆ New()

VtkImageDataToPointCloudFilter * VtkImageDataToPointCloudFilter::New ( )
static

Create a new objects (required because of VTKs reference counting)

References VtkImageDataToPointCloudFilter().

Referenced by VtkCompositeImageToPointCloudFilter::init(), and main().

◆ PrintSelf()

void VtkImageDataToPointCloudFilter::PrintSelf ( ostream & os,
vtkIndent indent )
override

Prints information about the filter.

Definition at line 25 of file VtkImageDataToPointCloudFilter.cpp.

26{
27 this->Superclass::PrintSelf(os, indent);
28}

◆ RequestData()

int VtkImageDataToPointCloudFilter::RequestData ( vtkInformation * request,
vtkInformationVector ** inputVector,
vtkInformationVector * outputVector )
overrideprotected

Updates the graphical object.

Definition at line 37 of file VtkImageDataToPointCloudFilter.cpp.

41{
42 vtkDebugMacro(<< "Executing VtkImageDataToPointCloudFilter");
43
44 vtkInformation* input_info = inputVector[0]->GetInformationObject(0);
45 vtkInformation* output_info = outputVector->GetInformationObject(0);
46 vtkImageData* input = vtkImageData::SafeDownCast(
47 input_info->Get(vtkDataObject::DATA_OBJECT()));
48 vtkPolyData* output = vtkPolyData::SafeDownCast(
49 output_info->Get(vtkDataObject::DATA_OBJECT()));
50
51 if (input->GetScalarSize() != sizeof(float))
52 {
53 vtkDebugMacro(
54 "Existing data does not have float-type and cannot be processed. "
55 "Aborting filter process...");
56 return 0;
57 }
58 float* pixvals = static_cast<float*>(input->GetScalarPointer());
59 int const n_comp = input->GetNumberOfScalarComponents();
60 if (n_comp < 1)
61 {
62 vtkDebugMacro("Error reading number of components.");
63 }
64 if (n_comp > 2)
65 {
66 vtkDebugMacro(
67 "RGB colours detected. Using only first channel for computation.");
68 }
69
70 vtkIdType const n_points = input->GetNumberOfPoints();
71 if (n_points == 0)
72 {
73 vtkDebugMacro("No data found!");
74 return -1;
75 }
76
77 double spacing[3];
78 input->GetSpacing(spacing);
79
80 if (MinValueRange == -1 || MaxValueRange == -1)
81 {
82 double range[2];
83 vtkPointData* input_data = input->GetPointData();
84 input_data->GetArray(0)->GetRange(range);
85 MinValueRange = range[0];
86 MaxValueRange = range[1];
87 }
88
89 std::vector<vtkIdType> density;
90 density.reserve(n_points);
91 for (vtkIdType i = 0; i < n_points; ++i)
92 {
93 // Skip transparent pixels
94 if (n_comp == 2 || n_comp == 4)
95 {
96 if (pixvals[(i + 1) * n_comp - 1] < 0.00000001f)
97 {
98 density.push_back(0);
99 continue;
100 }
101 }
102 double const val(pixvals[i * n_comp]);
103 double const calc_gamma = (IsLinear) ? 1 : Gamma;
104 std::size_t const pnts_per_cell =
105 interpolate(MinValueRange, MaxValueRange, val, calc_gamma);
106 density.push_back(static_cast<vtkIdType>(
107 std::floor(pnts_per_cell * GetPointScaleFactor())));
108 }
109
110 // Allocate the space needed for output objects
111 auto const sum =
112 std::accumulate(density.begin(), density.end(), vtkIdType{0});
113 vtkSmartPointer<vtkPoints> new_points = vtkSmartPointer<vtkPoints>::New();
114 new_points->SetNumberOfPoints(sum);
115 vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
116 cells->Allocate(sum);
117 vtkSmartPointer<vtkDoubleArray> intensity =
118 vtkSmartPointer<vtkDoubleArray>::New();
119 intensity->Allocate(sum);
120 intensity->SetName("Intensity");
121 double const half_cellsize(spacing[0] / 2.0);
122 std::size_t pnt_idx(0);
123 for (std::size_t i = 0; i < static_cast<std::size_t>(n_points); ++i)
124 {
125 if (density[i] == 0)
126 {
127 continue;
128 }
129 double p[3];
130 input->GetPoint(i, p);
131
132 MathLib::Point3d min_pnt{std::array<double, 3>{
133 {p[0] - half_cellsize, p[1] - half_cellsize, MinHeight}}};
134 MathLib::Point3d max_pnt{std::array<double, 3>{
135 {p[0] + half_cellsize, p[1] + half_cellsize, MaxHeight}}};
136 createPoints(new_points, cells, pnt_idx, density[i], min_pnt, max_pnt);
137 for (vtkIdType j = 0; j < density[i]; ++j)
138 {
139 intensity->InsertValue(pnt_idx + j, pixvals[i * n_comp]);
140 }
141 pnt_idx += density[i];
142 }
143
144 output->SetPoints(new_points);
145 output->SetVerts(cells);
146 output->GetPointData()->AddArray(intensity);
147 output->GetPointData()->SetActiveAttribute("Intensity",
148 vtkDataSetAttributes::SCALARS);
149 output->Squeeze();
150
151 vtkDebugMacro(<< "Created " << new_points->GetNumberOfPoints()
152 << " points.");
153 return 1;
154}
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.
std::size_t interpolate(double low, double high, double p, double gamma) const

References createPoints(), Gamma, interpolate(), IsLinear, MaxHeight, MaxValueRange, MinHeight, and MinValueRange.

◆ useLinearInterpolation()

void VtkImageDataToPointCloudFilter::useLinearInterpolation ( )
inline

Definition at line 27 of file VtkImageDataToPointCloudFilter.h.

27{ SetIsLinear(true); }

◆ useLogarithmicInterpolation()

void VtkImageDataToPointCloudFilter::useLogarithmicInterpolation ( double gamma)
inline

Definition at line 28 of file VtkImageDataToPointCloudFilter.h.

29 {
30 SetIsLinear(false);
31 SetGamma(gamma);
32 }

◆ vtkGetMacro() [1/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( Gamma ,
double  )

References Gamma.

◆ vtkGetMacro() [2/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( IsLinear ,
bool  )

References IsLinear.

◆ vtkGetMacro() [3/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MaxHeight ,
double  )

References MaxHeight.

◆ vtkGetMacro() [4/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MaxNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkGetMacro() [5/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MaxValueRange ,
double  )

References MaxValueRange.

◆ vtkGetMacro() [6/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MinHeight ,
double  )

References MinHeight.

◆ vtkGetMacro() [7/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MinNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkGetMacro() [8/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MinValueRange ,
double  )

References MinValueRange.

◆ vtkGetMacro() [9/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( PointScaleFactor ,
double  )

References PointScaleFactor.

◆ vtkSetMacro() [1/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( Gamma ,
double  )

References Gamma.

◆ vtkSetMacro() [2/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( IsLinear ,
bool  )

References IsLinear.

◆ vtkSetMacro() [3/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MaxHeight ,
double  )

References MaxHeight.

◆ vtkSetMacro() [4/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MaxNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkSetMacro() [5/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MaxValueRange ,
double  )

References MaxValueRange.

◆ vtkSetMacro() [6/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MinHeight ,
double  )

References MinHeight.

◆ vtkSetMacro() [7/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MinNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkSetMacro() [8/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MinValueRange ,
double  )

References MinValueRange.

◆ vtkSetMacro() [9/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( PointScaleFactor ,
double  )

References PointScaleFactor.

◆ vtkTypeMacro()

VtkImageDataToPointCloudFilter::vtkTypeMacro ( VtkImageDataToPointCloudFilter ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ Gamma

double VtkImageDataToPointCloudFilter::Gamma {1.0}
protected

Definition at line 72 of file VtkImageDataToPointCloudFilter.h.

72{1.0};

Referenced by RequestData(), vtkGetMacro(), and vtkSetMacro().

◆ IsLinear

bool VtkImageDataToPointCloudFilter::IsLinear {true}
protected

Definition at line 80 of file VtkImageDataToPointCloudFilter.h.

80{true};

Referenced by RequestData(), vtkGetMacro(), and vtkSetMacro().

◆ MaxHeight

double VtkImageDataToPointCloudFilter::MaxHeight {1000}
protected

Definition at line 77 of file VtkImageDataToPointCloudFilter.h.

77{1000};

Referenced by RequestData(), vtkGetMacro(), and vtkSetMacro().

◆ MaxNumberOfPointsPerCell

vtkIdType VtkImageDataToPointCloudFilter::MaxNumberOfPointsPerCell {20}
protected

Definition at line 79 of file VtkImageDataToPointCloudFilter.h.

79{20};

Referenced by interpolate(), vtkGetMacro(), and vtkSetMacro().

◆ MaxValueRange

double VtkImageDataToPointCloudFilter::MaxValueRange {-1}
protected

Definition at line 75 of file VtkImageDataToPointCloudFilter.h.

75{-1};

Referenced by RequestData(), vtkGetMacro(), and vtkSetMacro().

◆ MinHeight

double VtkImageDataToPointCloudFilter::MinHeight {0}
protected

Definition at line 76 of file VtkImageDataToPointCloudFilter.h.

76{0};

Referenced by RequestData(), vtkGetMacro(), and vtkSetMacro().

◆ MinNumberOfPointsPerCell

vtkIdType VtkImageDataToPointCloudFilter::MinNumberOfPointsPerCell {1}
protected

Definition at line 78 of file VtkImageDataToPointCloudFilter.h.

78{1};

Referenced by interpolate(), vtkGetMacro(), and vtkSetMacro().

◆ MinValueRange

double VtkImageDataToPointCloudFilter::MinValueRange {-1}
protected

Definition at line 74 of file VtkImageDataToPointCloudFilter.h.

74{-1};

Referenced by RequestData(), vtkGetMacro(), and vtkSetMacro().

◆ PointScaleFactor

double VtkImageDataToPointCloudFilter::PointScaleFactor {1.0}
protected

Definition at line 73 of file VtkImageDataToPointCloudFilter.h.

73{1.0};

Referenced by vtkGetMacro(), and vtkSetMacro().


The documentation for this class was generated from the following files: