VTK  9.3.0
vtkContinuousScatterplot.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
160 #ifndef vtkContinuousScatterplot_h
161 #define vtkContinuousScatterplot_h
162 
163 #include "vtkImageAlgorithm.h"
164 #include "vtkInfovisCoreModule.h" // For export macro
165 
166 VTK_ABI_NAMESPACE_BEGIN
167 class VTKINFOVISCORE_EXPORT vtkContinuousScatterplot : public vtkImageAlgorithm
168 {
169 public:
172  void PrintSelf(ostream& os, vtkIndent indent) override;
173 
177  vtkGetMacro(Epsilon, double);
178 
182  vtkSetMacro(Epsilon, double);
183 
188  void SetField1(const char* fieldName, vtkIdType ResX);
189 
194  void SetField2(const char* fieldName, vtkIdType ResY);
195 
196 protected:
198 
199  // Configure input port to accept only vtkUnstructuredGrid.
201 
202  // Configure out port to be a vtkImageData data set.
205 
206  // Set the tolerance used when comparing floating numbers for equality.
207  double Epsilon;
208 
209  // Names of the scalar fields to be used in the filter.
210  const char* Fields[2];
211 
212  // Resolution of the output image.
214 
215 private:
217  void operator=(const vtkContinuousScatterplot&) = delete;
218 };
219 VTK_ABI_NAMESPACE_END
220 #endif
Given a 3D domain space represented by an unstructured grid composed of tetrahedral cells with bivari...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int FillOutputPortInformation(int port, vtkInformation *info) override
These method should be reimplemented by subclasses that have more than a single input or single outpu...
void SetField1(const char *fieldName, vtkIdType ResX)
Specify the name of the first field to be used in subdividing the dataset.
void SetField2(const char *fieldName, vtkIdType ResY)
Specify the name of the second field to be used in subdividing the dataset.
static vtkContinuousScatterplot * New()
int FillInputPortInformation(int port, vtkInformation *info) override
These method should be reimplemented by subclasses that have more than a single input or single outpu...
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called in response to a REQUEST_DATA request from the executive.
Generic algorithm superclass for image algs.
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
int vtkIdType
Definition: vtkType.h:315