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 23 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

◆ ~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 163 of file VtkImageDataToPointCloudFilter.cpp.

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}
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 37 of file VtkImageDataToPointCloudFilter.cpp.

39{
40 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkImageData");
41 return 1;
42}

◆ getRandomNumber()

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

Returns a random number in [min, max].

Definition at line 183 of file VtkImageDataToPointCloudFilter.cpp.

185{
186 return (static_cast<double>(std::rand()) / RAND_MAX) * (max - min) + min;
187}

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 189 of file VtkImageDataToPointCloudFilter.cpp.

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}

References MaxNumberOfPointsPerCell, and MinNumberOfPointsPerCell.

Referenced by RequestData().

◆ New()

static VtkImageDataToPointCloudFilter * VtkImageDataToPointCloudFilter::New ( )
static

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

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

◆ PrintSelf()

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

Prints information about the filter.

Definition at line 32 of file VtkImageDataToPointCloudFilter.cpp.

33{
34 this->Superclass::PrintSelf(os, indent);
35}

◆ RequestData()

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

Updates the graphical object.

Definition at line 44 of file VtkImageDataToPointCloudFilter.cpp.

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}
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 34 of file VtkImageDataToPointCloudFilter.h.

34{ SetIsLinear(true); }

◆ useLogarithmicInterpolation()

void VtkImageDataToPointCloudFilter::useLogarithmicInterpolation ( double gamma)
inline

Definition at line 35 of file VtkImageDataToPointCloudFilter.h.

36 {
37 SetIsLinear(false);
38 SetGamma(gamma);
39 }

◆ vtkGetMacro() [1/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( Gamma ,
double  )

◆ vtkGetMacro() [2/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( IsLinear ,
bool  )

◆ vtkGetMacro() [3/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MaxHeight ,
double  )

◆ vtkGetMacro() [4/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MaxNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkGetMacro() [5/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MaxValueRange ,
double  )

◆ vtkGetMacro() [6/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MinHeight ,
double  )

◆ vtkGetMacro() [7/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MinNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkGetMacro() [8/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( MinValueRange ,
double  )

◆ vtkGetMacro() [9/9]

VtkImageDataToPointCloudFilter::vtkGetMacro ( PointScaleFactor ,
double  )

◆ vtkSetMacro() [1/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( Gamma ,
double  )

◆ vtkSetMacro() [2/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( IsLinear ,
bool  )

◆ vtkSetMacro() [3/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MaxHeight ,
double  )

◆ vtkSetMacro() [4/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MaxNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkSetMacro() [5/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MaxValueRange ,
double  )

◆ vtkSetMacro() [6/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MinHeight ,
double  )

◆ vtkSetMacro() [7/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MinNumberOfPointsPerCell ,
vtkIdType  )

◆ vtkSetMacro() [8/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( MinValueRange ,
double  )

◆ vtkSetMacro() [9/9]

VtkImageDataToPointCloudFilter::vtkSetMacro ( PointScaleFactor ,
double  )

◆ vtkTypeMacro()

VtkImageDataToPointCloudFilter::vtkTypeMacro ( VtkImageDataToPointCloudFilter ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ Gamma

double VtkImageDataToPointCloudFilter::Gamma {1.0}
protected

Definition at line 79 of file VtkImageDataToPointCloudFilter.h.

79{1.0};

Referenced by RequestData().

◆ IsLinear

bool VtkImageDataToPointCloudFilter::IsLinear {true}
protected

Definition at line 87 of file VtkImageDataToPointCloudFilter.h.

87{true};

Referenced by RequestData().

◆ MaxHeight

double VtkImageDataToPointCloudFilter::MaxHeight {1000}
protected

Definition at line 84 of file VtkImageDataToPointCloudFilter.h.

84{1000};

Referenced by RequestData().

◆ MaxNumberOfPointsPerCell

vtkIdType VtkImageDataToPointCloudFilter::MaxNumberOfPointsPerCell {20}
protected

Definition at line 86 of file VtkImageDataToPointCloudFilter.h.

86{20};

Referenced by interpolate().

◆ MaxValueRange

double VtkImageDataToPointCloudFilter::MaxValueRange {-1}
protected

Definition at line 82 of file VtkImageDataToPointCloudFilter.h.

82{-1};

Referenced by RequestData().

◆ MinHeight

double VtkImageDataToPointCloudFilter::MinHeight {0}
protected

Definition at line 83 of file VtkImageDataToPointCloudFilter.h.

83{0};

Referenced by RequestData().

◆ MinNumberOfPointsPerCell

vtkIdType VtkImageDataToPointCloudFilter::MinNumberOfPointsPerCell {1}
protected

Definition at line 85 of file VtkImageDataToPointCloudFilter.h.

85{1};

Referenced by interpolate().

◆ MinValueRange

double VtkImageDataToPointCloudFilter::MinValueRange {-1}
protected

Definition at line 81 of file VtkImageDataToPointCloudFilter.h.

81{-1};

Referenced by RequestData().

◆ PointScaleFactor

double VtkImageDataToPointCloudFilter::PointScaleFactor {1.0}
protected

Definition at line 80 of file VtkImageDataToPointCloudFilter.h.

80{1.0};

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