VTK  9.3.0
vtkPlotHistogram2D.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
3 
12 #ifndef vtkPlotHistogram2D_h
13 #define vtkPlotHistogram2D_h
14 
15 #include "vtkChartsCoreModule.h" // For export macro
16 #include "vtkPlot.h"
17 #include "vtkRect.h" // Needed for vtkRectf
18 #include "vtkSmartPointer.h" // Needed for SP ivars
19 
20 #include <string> // Needed for std::string
21 
22 VTK_ABI_NAMESPACE_BEGIN
23 
24 class vtkDataArray;
25 class vtkImageData;
26 class vtkScalarsToColors;
27 
28 class VTKCHARTSCORE_EXPORT vtkPlotHistogram2D : public vtkPlot
29 {
30 public:
31  vtkTypeMacro(vtkPlotHistogram2D, vtkPlot);
32  void PrintSelf(ostream& os, vtkIndent indent) override;
33 
38 
44  void Update() override;
45 
49  bool Paint(vtkContext2D* painter) override;
50 
55  virtual void SetInputData(vtkImageData* data, vtkIdType z = 0);
56  void SetInputData(vtkTable*) override {}
57  void SetInputData(vtkTable*, const vtkStdString&, const vtkStdString&) override {}
58 
63 
69 
74 
75  void GetBounds(double bounds[4]) override;
76 
77  virtual void SetPosition(const vtkRectf& pos);
78  virtual vtkRectf GetPosition();
79 
81 
86  vtkSetMacro(ArrayName, std::string);
87  vtkGetMacro(ArrayName, std::string);
89 
107  const vtkVector2d& plotPos, vtkIdType seriesIndex, vtkIdType segmentIndex) override;
108 
118  vtkVector2f* location, vtkIdType* segmentId) override;
120 
127  bool UpdateCache() override;
128 
129 protected:
132 
137 
138 private:
139  vtkPlotHistogram2D(const vtkPlotHistogram2D&) = delete;
140  void operator=(const vtkPlotHistogram2D&) = delete;
141 
146  static vtkIdType GetLabelIndexFromValue(double value, vtkAxis* axis);
151  static inline bool CanComputeMagnitude(int nbComponents);
152 
157  inline vtkDataArray* GetSelectedArray();
165  void* GetInputArrayPointer(int& nbComponents);
172  double GetInputArrayValue(int x, int y, int z);
173 
174  std::string ArrayName;
175  vtkSmartPointer<vtkDataArray> MagnitudeArray;
176 };
177 
178 VTK_ABI_NAMESPACE_END
179 #endif // vtkPlotHistogram2D_h
takes care of drawing 2D axes
Definition: vtkAxis.h:61
Class for drawing 2D primitives to a graphical context.
Definition: vtkContext2D.h:50
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
topologically and geometrically regular array of data
Definition: vtkImageData.h:43
a simple class to control print indentation
Definition: vtkIndent.h:29
vtkIdType GetNearestPoint(const vtkVector2f &point, const vtkVector2f &tolerance, vtkVector2f *location, vtkIdType *segmentId) override
Function to query a plot for the nearest point to the specified coordinate.
vtkScalarsToColors * GetTransferFunction()
Get the color transfer function that is used to generate the histogram.
bool Paint(vtkContext2D *painter) override
Paint event for the item, called whenever it needs to be drawn.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSmartPointer< vtkImageData > Output
void Update() override
Perform any updates to the item that may be necessary before rendering.
virtual vtkRectf GetPosition()
~vtkPlotHistogram2D() override
bool UpdateCache() override
Update the internal cache.
void SetInputData(vtkTable *, const vtkStdString &, const vtkStdString &) override
This is a convenience function to set the input table and the x, y column for the plot.
vtkStdString GetTooltipLabel(const vtkVector2d &plotPos, vtkIdType seriesIndex, vtkIdType segmentIndex) override
Generate and return the tooltip label string for this plot The segmentIndex parameter is ignored.
virtual void SetInputData(vtkImageData *data, vtkIdType z=0)
Set the input.
vtkImageData * GetInputImageData()
Get the input table used by the plot.
void SetInputData(vtkTable *) override
This is a convenience function to set the input table and the x, y column for the plot.
static vtkPlotHistogram2D * New()
Creates a new object.
void SetTransferFunction(vtkScalarsToColors *transfer)
Set the color transfer function that will be used to generate the 2D histogram.
virtual void SetPosition(const vtkRectf &pos)
vtkSmartPointer< vtkImageData > Input
void GetBounds(double bounds[4]) override
Get the bounds for this plot as (Xmin, Xmax, Ymin, Ymax).
vtkSmartPointer< vtkScalarsToColors > TransferFunction
Abstract class for 2D plots.
Definition: vtkPlot.h:44
virtual vtkIdType GetNearestPoint(const vtkVector2f &point, const vtkVector2f &tolerance, vtkVector2f *location, vtkIdType *segmentId)
Function to query a plot for the nearest point to the specified coordinate.
Superclass for mapping scalar values to colors.
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:29
A table, which contains similar-typed columns of data.
Definition: vtkTable.h:59
@ point
Definition: vtkX3D.h:236
@ location
Definition: vtkX3D.h:406
@ value
Definition: vtkX3D.h:220
@ data
Definition: vtkX3D.h:315
@ string
Definition: vtkX3D.h:490
int vtkIdType
Definition: vtkType.h:315