VTK  9.6.20260305
vtkLabeledImagePointSampler.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
72
73#ifndef vtkLabeledImagePointSampler_h
74#define vtkLabeledImagePointSampler_h
75
76#include "vtkContourValues.h" // Manage countour values
77#include "vtkFiltersMeshingModule.h" // For export macro
79
80VTK_ABI_NAMESPACE_BEGIN
81class vtkImageData;
82
83class VTKFILTERSMESHING_EXPORT vtkLabeledImagePointSampler : public vtkPolyDataAlgorithm
84{
85public:
87
93 void PrintSelf(ostream& os, vtkIndent indent) override;
95
97
111 vtkSetMacro(BackgroundPointMapping, int);
112 vtkGetMacro(BackgroundPointMapping, int);
113 vtkBooleanMacro(BackgroundPointMapping, int);
114 vtkSetMacro(BackgroundPointLabel, int);
115 vtkGetMacro(BackgroundPointLabel, int);
117
119
128
132 vtkSetClampMacro(DensityDistribution, int, LINEAR, EXPONENTIAL);
133 vtkGetMacro(DensityDistribution, int);
137
139
145 vtkSetClampMacro(N, unsigned int, 1, 100);
146 vtkGetMacro(N, unsigned int);
148
150
160 vtkSetMacro(Randomize, vtkTypeBool);
161 vtkGetMacro(Randomize, vtkTypeBool);
162 vtkBooleanMacro(Randomize, vtkTypeBool);
164
166
175 vtkSetVector2Macro(RandomProbabilityRange, double);
176 vtkGetVectorMacro(RandomProbabilityRange, double, 2);
178
180
189
195 vtkSetClampMacro(OutputType, int, ALL_POINTS, BACKGROUND_POINTS);
196 vtkGetMacro(OutputType, int);
201
203
212 vtkSetMacro(Joggle, vtkTypeBool);
213 vtkGetMacro(Joggle, vtkTypeBool);
214 vtkBooleanMacro(Joggle, vtkTypeBool);
215 vtkSetClampMacro(JoggleRadius, double, 0, 0.5);
216 vtkGetMacro(JoggleRadius, double);
221 vtkBooleanMacro(ConstrainJoggle, vtkTypeBool);
223 vtkSetClampMacro(JoggleConstraint, double, 0, 1);
224 vtkGetMacro(JoggleConstraint, double);
226
227 //---------------The following are methods used to set/get and generate
228 //---------------label values.
230
242 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
243 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
245
247
250 double GetValue(int i) { return this->Labels->GetValue(i); }
251 double GetLabel(int i) { return this->Labels->GetValue(i); }
253
255
259 double* GetValues() { return this->Labels->GetValues(); }
260 double* GetLabels() { return this->Labels->GetValues(); }
262
264
269 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
270 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
272
274
281 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
282 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
284
286
289 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
290 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
292
294
298 void GenerateLabels(int numLabels, double range[2])
299 {
300 this->Labels->GenerateValues(numLabels, range);
301 }
302 void GenerateValues(int numContours, double range[2])
303 {
304 this->Labels->GenerateValues(numContours, range);
305 }
306 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
307 {
308 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
309 }
310 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
311 {
312 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
313 }
314
315 //---------------Done defining label-related methods.
316
318
326 vtkBooleanMacro(GenerateVerts, vtkTypeBool);
328
333
334protected:
336 ~vtkLabeledImagePointSampler() override = default;
337
341 unsigned int N;
343
346
352
354
356
358 int FillInputPortInformation(int port, vtkInformation* info) override;
359
360private:
362 void operator=(const vtkLabeledImagePointSampler&) = delete;
363};
364
365VTK_ABI_NAMESPACE_END
366#endif
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
~vtkLabeledImagePointSampler() override=default
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void SetNumberOfContours(int number)
Set the number of labels to place into the list.
double * GetLabels()
Get a pointer to an array of labels.
double GetValue(int i)
Get the ith label value.
void GenerateValues(int numContours, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
DensityDistributionType
Used to indicate the type of point density distribution.
void GetValues(double *contourValues)
Fill a supplied list with label values.
void SetOutputTypeToAllPoints()
Used to indicate which generated points to output.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, obtaining type information, and printing information.
void SetDensityDistributionToLinear()
Used to indicate the type of point density distribution.
OutputLabelsSelection
Used to indicate which generated points to output.
void SetOutputTypeToLabeledPoints()
Used to indicate which generated points to output.
double GetLabel(int i)
Get the ith label value.
virtual void SetOutputType(int)
Indicate what points to output.
virtual void SetDensityDistribution(int)
Specify the type of point selection distribution (i.e., the specification of the point density).
vtkSmartPointer< vtkContourValues > Labels
double * GetValues()
Get a pointer to an array of labels.
vtkMTimeType GetMTime() override
Because we delegate to vtkContourValues.
void SetDensityDistributionToExponential()
Used to indicate the type of point density distribution.
void SetValue(int i, double value)
Set a particular label value at label number i.
void SetLabel(int i, double value)
Set a particular label value at label number i.
static vtkLabeledImagePointSampler * New()
Standard methods for instantiation, obtaining type information, and printing information.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void SetOutputTypeToBackgroundPoints()
Used to indicate which generated points to output.
vtkIdType GetNumberOfLabels()
Get the number of labels in the list of label values.
void GetLabels(double *contourValues)
Fill a supplied list with label values.
void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
void GenerateLabels(int numLabels, double range[2])
Generate numLabels equally spaced labels between the specified range.
void GenerateValues(int numContours, double range[2])
Generate numLabels equally spaced labels between the specified range.
void SetNumberOfLabels(int number)
Set the number of labels to place into the list.
vtkIdType GetNumberOfContours()
Get the number of labels in the list of label values.
Hold a reference to a vtkObjectBase instance.
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:363
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318