VTK  9.3.0
vtkPointInterpolator.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
65#ifndef vtkPointInterpolator_h
66#define vtkPointInterpolator_h
67
68#include "vtkDataSetAlgorithm.h"
69#include "vtkFiltersPointsModule.h" // For export macro
70#include "vtkStdString.h" // For vtkStdString ivars
71#include <vector> //For STL vector
72
73VTK_ABI_NAMESPACE_BEGIN
75class vtkIdList;
76class vtkDoubleArray;
78class vtkCharArray;
79
80class VTKFILTERSPOINTS_EXPORT vtkPointInterpolator : public vtkDataSetAlgorithm
81{
82public:
84
90 void PrintSelf(ostream& os, vtkIndent indent) override;
92
94
104
112
114
120 vtkGetObjectMacro(Locator, vtkAbstractPointLocator);
122
124
130 vtkGetObjectMacro(Kernel, vtkInterpolationKernel);
132
134 {
135 MASK_POINTS = 0,
136 NULL_VALUE = 1,
137 CLOSEST_POINT = 2
138 };
139
141
152 vtkSetMacro(NullPointsStrategy, int);
153 vtkGetMacro(NullPointsStrategy, int);
154 void SetNullPointsStrategyToMaskPoints() { this->SetNullPointsStrategy(MASK_POINTS); }
155 void SetNullPointsStrategyToNullValue() { this->SetNullPointsStrategy(NULL_VALUE); }
156 void SetNullPointsStrategyToClosestPoint() { this->SetNullPointsStrategy(CLOSEST_POINT); }
158
160
166 vtkSetMacro(ValidPointsMaskArrayName, vtkStdString);
167 vtkGetMacro(ValidPointsMaskArrayName, vtkStdString);
169
171
176 vtkSetMacro(NullValue, double);
177 vtkGetMacro(NullValue, double);
179
181
185 void AddExcludedArray(const vtkStdString& excludedArray)
186 {
187 this->ExcludedArrays.push_back(excludedArray);
188 this->Modified();
189 }
191
193
197 {
198 this->ExcludedArrays.clear();
199 this->Modified();
200 }
202
206 int GetNumberOfExcludedArrays() { return static_cast<int>(this->ExcludedArrays.size()); }
207
209
212 const char* GetExcludedArray(int i)
213 {
214 if (i < 0 || i >= static_cast<int>(this->ExcludedArrays.size()))
215 {
216 return nullptr;
217 }
218 return this->ExcludedArrays[i].c_str();
219 }
221
223
229 vtkSetMacro(PromoteOutputArrays, bool);
230 vtkBooleanMacro(PromoteOutputArrays, bool);
231 vtkGetMacro(PromoteOutputArrays, bool);
233
235
239 vtkSetMacro(PassPointArrays, bool);
240 vtkBooleanMacro(PassPointArrays, bool);
241 vtkGetMacro(PassPointArrays, bool);
243
245
249 vtkSetMacro(PassCellArrays, bool);
250 vtkBooleanMacro(PassCellArrays, bool);
251 vtkGetMacro(PassCellArrays, bool);
253
255
259 vtkSetMacro(PassFieldArrays, bool);
260 vtkBooleanMacro(PassFieldArrays, bool);
261 vtkGetMacro(PassFieldArrays, bool);
263
268
269protected:
272
275
277 double NullValue;
280
281 std::vector<vtkStdString> ExcludedArrays;
282
284
288
292
296 virtual void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
297
303
308 vtkImageData* input, int dims[3], double origin[3], double spacing[3]);
309
310private:
312 void operator=(const vtkPointInterpolator&) = delete;
313};
314
315VTK_ABI_NAMESPACE_END
316#endif
abstract class to quickly locate points in 3-space
Proxy object to connect input/output ports.
dynamic, self-adjusting array of char
general representation of visualization data
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition vtkDataSet.h:62
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:32
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
base class for interpolation kernels
virtual void Modified()
Update the modification time for this object.
interpolate over point cloud using various kernels
vtkInterpolationKernel * Kernel
int GetNumberOfExcludedArrays()
Return the number of excluded arrays.
std::vector< vtkStdString > ExcludedArrays
void SetNullPointsStrategyToNullValue()
Specify a strategy to use when encountering a "null" point during the interpolation process.
void AddExcludedArray(const vtkStdString &excludedArray)
Adds an array to the list of arrays which are to be excluded from the interpolation process.
vtkStdString ValidPointsMaskArrayName
void SetSourceData(vtkDataObject *source)
Specify the dataset Pc that will be probed by the input points P.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiating, obtaining type information, and printing.
vtkDataObject * GetSource()
Specify the dataset Pc that will be probed by the input points P.
virtual void PassAttributeData(vtkDataSet *input, vtkDataObject *source, vtkDataSet *output)
Call at end of RequestData() to pass attribute data respecting the PassCellArrays,...
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the dataset Pc that will be probed by the input points P.
static vtkPointInterpolator * New()
Standard methods for instantiating, obtaining type information, and printing.
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when each filter in the pipeline decides what portion of its inp...
virtual void Probe(vtkDataSet *input, vtkDataSet *source, vtkDataSet *output)
Virtual for specialized subclass(es)
void SetNullPointsStrategyToClosestPoint()
Specify a strategy to use when encountering a "null" point during the interpolation process.
vtkAbstractPointLocator * Locator
void SetKernel(vtkInterpolationKernel *kernel)
Specify an interpolation kernel.
void ExtractImageDescription(vtkImageData *input, int dims[3], double origin[3], double spacing[3])
Internal method to extract image metadata.
~vtkPointInterpolator() override
void ClearExcludedArrays()
Clears the contents of excluded array list.
void SetNullPointsStrategyToMaskPoints()
Specify a strategy to use when encountering a "null" point during the interpolation process.
const char * GetExcludedArray(int i)
Return the name of the ith excluded array.
vtkMTimeType GetMTime() override
Get the MTime of this object also considering the locator and kernel.
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks for Information.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
void SetLocator(vtkAbstractPointLocator *locator)
Specify a point locator.
Wrapper around std::string to keep symbols short.
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270