OGS
VtkImageDataToSurfacePointsFilter Class Reference

Detailed Description

A VTK filter that creates a point cloud representing a surface defined by pixel values

Definition at line 23 of file VtkImageDataToSurfacePointsFilter.h.

#include <VtkImageDataToSurfacePointsFilter.h>

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

Public Member Functions

 vtkTypeMacro (VtkImageDataToSurfacePointsFilter, vtkPolyDataAlgorithm)
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Prints information about itself.
 
 vtkGetMacro (PointsPerPixel, vtkIdType)
 
 vtkSetMacro (PointsPerPixel, vtkIdType)
 

Static Public Member Functions

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

Protected Member Functions

 VtkImageDataToSurfacePointsFilter ()
 
 ~VtkImageDataToSurfacePointsFilter () 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.
 

Private Member Functions

void createPointSurface (vtkSmartPointer< vtkPoints > &points, vtkSmartPointer< vtkCellArray > &cells, std::size_t pnt_idx, MathLib::Point3d const &min_pnt, MathLib::Point3d const &max_pnt, GeoLib::Raster const &raster)
 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].
 

Private Attributes

vtkIdType PointsPerPixel {20}
 

Constructor & Destructor Documentation

◆ VtkImageDataToSurfacePointsFilter()

VtkImageDataToSurfacePointsFilter::VtkImageDataToSurfacePointsFilter ( )
protecteddefault

◆ ~VtkImageDataToSurfacePointsFilter()

VtkImageDataToSurfacePointsFilter::~VtkImageDataToSurfacePointsFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ createPointSurface()

void VtkImageDataToSurfacePointsFilter::createPointSurface ( vtkSmartPointer< vtkPoints > & points,
vtkSmartPointer< vtkCellArray > & cells,
std::size_t pnt_idx,
MathLib::Point3d const & min_pnt,
MathLib::Point3d const & max_pnt,
GeoLib::Raster const & raster )
private

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

Definition at line 146 of file VtkImageDataToSurfacePointsFilter.cpp.

153{
154 auto const n_points(static_cast<std::size_t>(this->GetPointsPerPixel()));
155 for (std::size_t i = 0; i < n_points; ++i)
156 {
157 double p[3];
158 p[0] = getRandomNumber(min_pnt[0], max_pnt[0]);
159 p[1] = getRandomNumber(min_pnt[1], max_pnt[1]);
160 p[2] = raster.interpolateValueAtPoint(
161 MathLib::Point3d(std::array<double, 3>{{p[0], p[1], 0}}));
162 points->SetPoint(pnt_idx + i, p);
163 cells->InsertNextCell(1);
164 cells->InsertCellPoint(pnt_idx + i);
165 }
166}
double getRandomNumber(double const &min, double const &max) const
Returns a random number in [min, max].

References getRandomNumber(), and GeoLib::Raster::interpolateValueAtPoint().

Referenced by RequestData().

◆ FillInputPortInformation()

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

Sets input port to vtkImageData.

Definition at line 37 of file VtkImageDataToSurfacePointsFilter.cpp.

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

◆ getRandomNumber()

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

Returns a random number in [min, max].

Definition at line 168 of file VtkImageDataToSurfacePointsFilter.cpp.

170{
171 return (static_cast<double>(std::rand()) / RAND_MAX) * (max - min) + min;
172}

Referenced by createPointSurface().

◆ New()

static VtkImageDataToSurfacePointsFilter * VtkImageDataToSurfacePointsFilter::New ( )
static

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

Referenced by VtkCompositeImageToSurfacePointsFilter::init().

◆ PrintSelf()

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

Prints information about itself.

Definition at line 32 of file VtkImageDataToSurfacePointsFilter.cpp.

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

◆ RequestData()

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

Updates the graphical object.

Definition at line 44 of file VtkImageDataToSurfacePointsFilter.cpp.

48{
49 vtkDebugMacro(<< "Executing VtkImageDataToSurfacePointsFilter");
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 void* pixvals = input->GetScalarPointer();
59 int 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 auto const n_points = static_cast<std::size_t>(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 int dimensions[3];
80 input->GetDimensions(dimensions);
81 double origin[3];
82 input->GetOrigin(origin);
83 MathLib::Point3d const ll(
84 std::array<double, 3>{{origin[0], origin[1], origin[2]}});
85
86 std::vector<double> pixels;
87 pixels.reserve(n_points);
88 for (std::size_t i = 0; i < n_points; ++i)
89 {
90 if ((n_comp == 2 || n_comp == 4) &&
91 ((static_cast<float*>(pixvals))[(i + 1) * n_comp - 1] <
92 0.00000001f))
93 {
94 pixels.push_back(-9999);
95 }
96 else
97 {
98 pixels.push_back((static_cast<float*>(pixvals))[i * n_comp]);
99 }
100 }
101 GeoLib::Raster const* const raster(new GeoLib::Raster(
102 {static_cast<std::size_t>(dimensions[0]),
103 static_cast<std::size_t>(dimensions[1]), 1, ll, spacing[0], -9999},
104 pixels.begin(),
105 pixels.end()));
106
107 vtkSmartPointer<vtkPoints> new_points = vtkSmartPointer<vtkPoints>::New();
108 new_points->SetNumberOfPoints(PointsPerPixel * n_points);
109 vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
110 cells->Allocate(PointsPerPixel * n_points);
111
112 double const half_cellsize(spacing[0] / 2.0);
113 std::size_t pnt_idx(0);
114 for (std::size_t i = 0; i < static_cast<std::size_t>(n_points); ++i)
115 {
116 // Skip transparent pixels
117 if (n_comp == 2 || n_comp == 4)
118 {
119 if ((static_cast<float*>(pixvals))[(i + 1) * n_comp - 1] <
120 0.00000001f)
121 {
122 continue;
123 }
124 }
125
126 double p[3];
127 input->GetPoint(i, p);
128 MathLib::Point3d min_pnt{std::array<double, 3>{
129 {p[0] - half_cellsize, p[1] - half_cellsize, 0}}};
130 MathLib::Point3d max_pnt{std::array<double, 3>{
131 {p[0] + half_cellsize, p[1] + half_cellsize, 0}}};
133 new_points, cells, pnt_idx, min_pnt, max_pnt, *raster);
134 pnt_idx += PointsPerPixel;
135 }
136
137 output->SetPoints(new_points);
138 output->SetVerts(cells);
139 output->Squeeze();
140
141 vtkDebugMacro(<< "Created " << new_points->GetNumberOfPoints()
142 << " points.");
143 return 1;
144}
Class Raster is used for managing raster data.
Definition Raster.h:49
void createPointSurface(vtkSmartPointer< vtkPoints > &points, vtkSmartPointer< vtkCellArray > &cells, std::size_t pnt_idx, MathLib::Point3d const &min_pnt, MathLib::Point3d const &max_pnt, GeoLib::Raster const &raster)
Creates the point objects based on the parameters set by the user.

References createPointSurface(), and PointsPerPixel.

◆ vtkGetMacro()

VtkImageDataToSurfacePointsFilter::vtkGetMacro ( PointsPerPixel ,
vtkIdType  )

◆ vtkSetMacro()

VtkImageDataToSurfacePointsFilter::vtkSetMacro ( PointsPerPixel ,
vtkIdType  )

◆ vtkTypeMacro()

VtkImageDataToSurfacePointsFilter::vtkTypeMacro ( VtkImageDataToSurfacePointsFilter ,
vtkPolyDataAlgorithm  )

Member Data Documentation

◆ PointsPerPixel

vtkIdType VtkImageDataToSurfacePointsFilter::PointsPerPixel {20}
private

Definition at line 60 of file VtkImageDataToSurfacePointsFilter.h.

60{20};

Referenced by RequestData().


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