VTK  9.3.20240419
vtkPairwiseExtractHistogram2D.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright 2009 Sandia Corporation
3 // SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-USGov
29 #ifndef vtkPairwiseExtractHistogram2D_h
30 #define vtkPairwiseExtractHistogram2D_h
31 
32 #include "vtkFiltersImagingModule.h" // For export macro
33 #include "vtkSmartPointer.h" //needed for smart pointer ivars
34 #include "vtkStatisticsAlgorithm.h"
35 VTK_ABI_NAMESPACE_BEGIN
36 class vtkCollection;
38 class vtkImageData;
39 class vtkIdTypeArray;
41 
42 class VTKFILTERSIMAGING_EXPORT vtkPairwiseExtractHistogram2D : public vtkStatisticsAlgorithm
43 {
44 public:
47  void PrintSelf(ostream& os, vtkIndent indent) override;
48 
50 
53  vtkSetVector2Macro(NumberOfBins, int);
54  vtkGetVector2Macro(NumberOfBins, int);
56 
58 
63  vtkSetMacro(CustomColumnRangeIndex, int);
64  void SetCustomColumnRangeByIndex(double, double);
66 
68 
73  void SetCustomColumnRange(int col, double range[2]);
74  void SetCustomColumnRange(int col, double rmin, double rmax);
76 
78 
81  vtkSetMacro(ScalarType, int);
82  void SetScalarTypeToUnsignedInt() { this->SetScalarType(VTK_UNSIGNED_INT); }
83  void SetScalarTypeToUnsignedLong() { this->SetScalarType(VTK_UNSIGNED_LONG); }
84  void SetScalarTypeToUnsignedShort() { this->SetScalarType(VTK_UNSIGNED_SHORT); }
85  void SetScalarTypeToUnsignedChar() { this->SetScalarType(VTK_UNSIGNED_CHAR); }
86  vtkGetMacro(ScalarType, int);
88 
92  double GetMaximumBinCount(int idx);
93 
98 
103  int GetBinRange(int idx, vtkIdType binX, vtkIdType binY, double range[4]);
104 
109  int GetBinRange(int idx, vtkIdType bin, double range[4]);
110 
115  void GetBinWidth(int idx, double bw[2]);
116 
121  double* GetHistogramExtents(int idx);
122 
127 
132 
134  {
135  HISTOGRAM_IMAGE = 3
136  };
137 
142 
143 protected:
146 
147  int NumberOfBins[2];
150 
153  class Internals;
154  Internals* Implementation;
155 
160  void Learn(vtkTable* inData, vtkTable* inParameters, vtkMultiBlockDataSet* outMeta) override;
161 
165  void Derive(vtkMultiBlockDataSet*) override {}
166 
171 
176 
180  void SelectAssessFunctor(vtkTable* vtkNotUsed(outData), vtkDataObject* vtkNotUsed(inMeta),
181  vtkStringArray* vtkNotUsed(rowNames), AssessFunctor*& vtkNotUsed(dfunc)) override
182  {
183  }
184 
189 
191 
193 
194 private:
196  void operator=(const vtkPairwiseExtractHistogram2D&) = delete;
197 };
198 
199 VTK_ABI_NAMESPACE_END
200 #endif
create and manipulate ordered lists of objects
Definition: vtkCollection.h:46
maintain an unordered list of data objects
general representation of visualization data
compute a 2D histogram between two columns of an input vtkTable.
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
Definition: vtkImageData.h:156
a simple class to control print indentation
Definition: vtkIndent.h:108
Store vtkAlgorithm input/output information.
Composite dataset that organizes datasets into blocks.
compute a 2D histogram between all adjacent columns of an input vtkTable.
void SetScalarTypeToUnsignedLong()
Set the scalar type for each of the computed histograms.
void Learn(vtkTable *inData, vtkTable *inParameters, vtkMultiBlockDataSet *outMeta) override
Execute the calculations required by the Learn option.
void SetCustomColumnRange(int col, double range[2])
More standard way to set the custom range for a particular column.
void SetCustomColumnRangeByIndex(double, double)
Strange method for setting an index to be used for setting custom column range.
void SetScalarTypeToUnsignedInt()
Set the scalar type for each of the computed histograms.
int GetBinRange(int idx, vtkIdType binX, vtkIdType binY, double range[4])
Compute the range of the bin located at position (binX,binY) in the 2D histogram at idx.
int GetBinRange(int idx, vtkIdType bin, double range[4])
Get the range of the of the bin located at 1D position index bin in the 2D histogram array at idx.
void Aggregate(vtkDataObjectCollection *, vtkMultiBlockDataSet *) override
Given a collection of models, calculate aggregate model.
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
void Derive(vtkMultiBlockDataSet *) override
Execute the calculations required by the Derive option.
void GetBinWidth(int idx, double bw[2])
Get the width of all of the bins.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSmartPointer< vtkIdTypeArray > OutputOutlierIds
~vtkPairwiseExtractHistogram2D() override
void SetScalarTypeToUnsignedShort()
Set the scalar type for each of the computed histograms.
void SetCustomColumnRange(int col, double rmin, double rmax)
More standard way to set the custom range for a particular column.
void SetScalarTypeToUnsignedChar()
Set the scalar type for each of the computed histograms.
static vtkPairwiseExtractHistogram2D * New()
double GetMaximumBinCount()
Get the maximum bin count over all histograms.
void Test(vtkTable *, vtkMultiBlockDataSet *, vtkTable *) override
Execute the calculations required by the Test option.
virtual vtkExtractHistogram2D * NewHistogramFilter()
Generate a new histogram filter.
vtkSmartPointer< vtkCollection > HistogramFilters
vtkExtractHistogram2D * GetHistogramFilter(int idx)
Get a pointer to the idx'th histogram filter.
vtkImageData * GetOutputHistogramImage(int idx)
Get the vtkImageData output of the idx'th histogram filter.
void SelectAssessFunctor(vtkTable *vtkNotUsed(outData), vtkDataObject *vtkNotUsed(inMeta), vtkStringArray *vtkNotUsed(rowNames), AssessFunctor *&vtkNotUsed(dfunc)) override
Provide the appropriate assessment functor.
double * GetHistogramExtents(int idx)
Get the histogram extents currently in use, either computed or set by the user for the idx'th histogr...
double GetMaximumBinCount(int idx)
Get the maximum bin count for a single histogram.
void Assess(vtkTable *, vtkMultiBlockDataSet *, vtkTable *) override
Execute the assess option.
A base class for a functor that assesses data.
Base class for statistics algorithms.
a vtkAbstractArray subclass for strings
A table, which contains similar-typed columns of data.
Definition: vtkTable.h:168
record modification and/or execution time
Definition: vtkTimeStamp.h:44
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
@ range
Definition: vtkX3D.h:238
int vtkIdType
Definition: vtkType.h:315
#define VTK_UNSIGNED_INT
Definition: vtkType.h:39
#define VTK_UNSIGNED_CHAR
Definition: vtkType.h:35
#define VTK_UNSIGNED_SHORT
Definition: vtkType.h:37
#define VTK_UNSIGNED_LONG
Definition: vtkType.h:41