VTK  9.1.0
vtkPairwiseExtractHistogram2D.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkPairwiseExtractHistogram2D.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
15/*-------------------------------------------------------------------------
16 Copyright 2009 Sandia Corporation.
17 Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
18 the U.S. Government retains certain rights in this software.
19-------------------------------------------------------------------------*/
45#ifndef vtkPairwiseExtractHistogram2D_h
46#define vtkPairwiseExtractHistogram2D_h
47
48#include "vtkFiltersImagingModule.h" // For export macro
49#include "vtkSmartPointer.h" //needed for smart pointer ivars
51class vtkCollection;
53class vtkImageData;
54class vtkIdTypeArray;
56
57class VTKFILTERSIMAGING_EXPORT vtkPairwiseExtractHistogram2D : public vtkStatisticsAlgorithm
58{
59public:
62 void PrintSelf(ostream& os, vtkIndent indent) override;
63
65
68 vtkSetVector2Macro(NumberOfBins, int);
69 vtkGetVector2Macro(NumberOfBins, int);
71
73
78 vtkSetMacro(CustomColumnRangeIndex, int);
79 void SetCustomColumnRangeByIndex(double, double);
81
83
88 void SetCustomColumnRange(int col, double range[2]);
89 void SetCustomColumnRange(int col, double rmin, double rmax);
91
93
96 vtkSetMacro(ScalarType, int);
97 void SetScalarTypeToUnsignedInt() { this->SetScalarType(VTK_UNSIGNED_INT); }
98 void SetScalarTypeToUnsignedLong() { this->SetScalarType(VTK_UNSIGNED_LONG); }
99 void SetScalarTypeToUnsignedShort() { this->SetScalarType(VTK_UNSIGNED_SHORT); }
100 void SetScalarTypeToUnsignedChar() { this->SetScalarType(VTK_UNSIGNED_CHAR); }
101 vtkGetMacro(ScalarType, int);
103
107 double GetMaximumBinCount(int idx);
108
113
118 int GetBinRange(int idx, vtkIdType binX, vtkIdType binY, double range[4]);
119
124 int GetBinRange(int idx, vtkIdType bin, double range[4]);
125
130 void GetBinWidth(int idx, double bw[2]);
131
136 double* GetHistogramExtents(int idx);
137
142
147
149 {
150 HISTOGRAM_IMAGE = 3
151 };
152
157
158protected:
161
162 int NumberOfBins[2];
165
168 class Internals;
169 Internals* Implementation;
170
175 void Learn(vtkTable* inData, vtkTable* inParameters, vtkMultiBlockDataSet* outMeta) override;
176
180 void Derive(vtkMultiBlockDataSet*) override {}
181
186
190 void Test(vtkTable*, vtkMultiBlockDataSet*, vtkTable*) override { return; }
191
195 void SelectAssessFunctor(vtkTable* vtkNotUsed(outData), vtkDataObject* vtkNotUsed(inMeta),
196 vtkStringArray* vtkNotUsed(rowNames), AssessFunctor*& vtkNotUsed(dfunc)) override
197 {
198 }
199
204
205 int FillOutputPortInformation(int port, vtkInformation* info) override;
206
208
209private:
211 void operator=(const vtkPairwiseExtractHistogram2D&) = delete;
212};
213
214#endif
create and manipulate ordered lists of objects
Definition: vtkCollection.h:53
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:157
a simple class to control print indentation
Definition: vtkIndent.h:113
Store vtkAlgorithm input/output information.
Composite dataset that organizes datasets into blocks.
compute a 2D histogram between all adjacent columns of an input vtkTable.
double * GetHistogramExtents(int idx)
Get the histogram extents currently in use, either computed or set by the user for the idx'th histogr...
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.
virtual vtkExtractHistogram2D * NewHistogramFilter()
Generate a new histogram filter.
void Aggregate(vtkDataObjectCollection *, vtkMultiBlockDataSet *) override
Given a collection of models, calculate aggregate model.
vtkImageData * GetOutputHistogramImage(int idx)
Get the vtkImageData output of the idx'th histogram filter.
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
vtkExtractHistogram2D * GetHistogramFilter(int idx)
Get a pointer to the idx'th histogram filter.
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.
vtkSmartPointer< vtkCollection > HistogramFilters
void SelectAssessFunctor(vtkTable *vtkNotUsed(outData), vtkDataObject *vtkNotUsed(inMeta), vtkStringArray *vtkNotUsed(rowNames), AssessFunctor *&vtkNotUsed(dfunc)) override
Provide the appropriate assessment functor.
double GetMaximumBinCount(int idx)
Get the maximum bin count for a single histogram.
void Assess(vtkTable *, vtkMultiBlockDataSet *, vtkTable *) override
Execute the assess option.
Hold a reference to a vtkObjectBase instance.
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:172
record modification and/or execution time
Definition: vtkTimeStamp.h:52
int vtkIdType
Definition: vtkType.h:332
#define VTK_UNSIGNED_INT
Definition: vtkType.h:51
#define VTK_UNSIGNED_CHAR
Definition: vtkType.h:47
#define VTK_UNSIGNED_SHORT
Definition: vtkType.h:49
#define VTK_UNSIGNED_LONG
Definition: vtkType.h:53