VTK  9.1.0
vtkGradientFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkGradientFilter.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 (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
46 #ifndef vtkGradientFilter_h
47 #define vtkGradientFilter_h
48 
49 #include "vtkDataSetAlgorithm.h"
50 #include "vtkFiltersGeneralModule.h" // For export macro
51 
52 class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
53 {
54 public:
56 
61  void PrintSelf(ostream& os, vtkIndent indent) override;
63 
66  {
67  All = 0,
68  Patch = 1,
69  DataSetMax = 2
70  };
71 
75  {
76  Zero = 0,
77  NaN = 1,
78  DataTypeMin = 2,
79  DataTypeMax = 3
80  };
81 
83 
89  virtual void SetInputScalars(int fieldAssociation, const char* name);
90  virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
92 
94 
99  vtkGetStringMacro(ResultArrayName);
100  vtkSetStringMacro(ResultArrayName);
102 
104 
109  vtkGetStringMacro(DivergenceArrayName);
110  vtkSetStringMacro(DivergenceArrayName);
112 
114 
119  vtkGetStringMacro(VorticityArrayName);
120  vtkSetStringMacro(VorticityArrayName);
122 
124 
129  vtkGetStringMacro(QCriterionArrayName);
130  vtkSetStringMacro(QCriterionArrayName);
132 
134 
143  vtkGetMacro(FasterApproximation, vtkTypeBool);
144  vtkSetMacro(FasterApproximation, vtkTypeBool);
145  vtkBooleanMacro(FasterApproximation, vtkTypeBool);
147 
149 
154  vtkSetMacro(ComputeGradient, vtkTypeBool);
155  vtkGetMacro(ComputeGradient, vtkTypeBool);
156  vtkBooleanMacro(ComputeGradient, vtkTypeBool);
158 
160 
166  vtkSetMacro(ComputeDivergence, vtkTypeBool);
167  vtkGetMacro(ComputeDivergence, vtkTypeBool);
168  vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
170 
172 
178  vtkSetMacro(ComputeVorticity, vtkTypeBool);
179  vtkGetMacro(ComputeVorticity, vtkTypeBool);
180  vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
182 
184 
191  vtkSetMacro(ComputeQCriterion, vtkTypeBool);
192  vtkGetMacro(ComputeQCriterion, vtkTypeBool);
193  vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
195 
197 
201  vtkSetClampMacro(ContributingCellOption, int, 0, 2);
202  vtkGetMacro(ContributingCellOption, int);
204 
206 
211  vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
212  vtkGetMacro(ReplacementValueOption, int);
214 
215 protected:
217  ~vtkGradientFilter() override;
218 
221 
227  virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
228  vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
229  vtkDataSet* output);
230 
236  virtual int ComputeRegularGridGradient(vtkDataArray* Array, int fieldAssociation,
237  bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output);
238 
246 
252 
258 
264 
270 
281 
287 
294 
301 
308 
314 
321 
322 private:
323  vtkGradientFilter(const vtkGradientFilter&) = delete;
324  void operator=(const vtkGradientFilter&) = delete;
325 };
326 
327 #endif //_vtkGradientFilter_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
A general filter for gradient estimation.
char * QCriterionArrayName
If non-null then it contains the name of the outputted Q criterion array.
~vtkGradientFilter() override
int ReplacementValueOption
Option to specify what replacement value or entities that don't have any gradient computed over them ...
char * VorticityArrayName
If non-null then it contains the name of the outputted vorticity array.
virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType)
These are basically a convenience method that calls SetInputArrayToProcess to set the array used as t...
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, obtaining type information, and printing.
vtkTypeBool ComputeVorticity
Flag to indicate that vorticity/curl of the input vector is to be computed.
vtkTypeBool ComputeGradient
Flag to indicate that the gradient of the input vector is to be computed.
int ContributingCellOption
Option to specify what cells to include in the gradient computation.
vtkTypeBool ComputeDivergence
Flag to indicate that the divergence of the input vector is to be computed.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
char * ResultArrayName
If non-null then it contains the name of the outputted gradient array.
ContributingCellEnum
Options to choose what cells contribute to the gradient calculation.
static vtkGradientFilter * New()
Standard methods for instantiation, obtaining type information, and printing.
int GetOutputArrayType(vtkDataArray *inputArray)
Get the proper array type to compute requested derivative quantities for.
virtual int ComputeUnstructuredGridGradient(vtkDataArray *Array, int fieldAssociation, vtkDataSet *input, bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet *output)
Compute the gradients for grids that are not a vtkImageData, vtkRectilinearGrid, or vtkStructuredGrid...
char * DivergenceArrayName
If non-null then it contains the name of the outputted divergence array.
virtual void SetInputScalars(int fieldAssociation, const char *name)
These are basically a convenience method that calls SetInputArrayToProcess to set the array used as t...
virtual int ComputeRegularGridGradient(vtkDataArray *Array, int fieldAssociation, bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet *output)
Compute the gradients for either a vtkImageData, vtkRectilinearGrid or a vtkStructuredGrid.
vtkTypeBool FasterApproximation
When this flag is on (default is off), the gradient filter will provide a less accurate (but close) a...
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when each filter in the pipeline decides what portion of its inp...
ReplacementValueEnum
The replacement value or entities that don't have any gradient computed over them based on the Contri...
vtkTypeBool ComputeQCriterion
Flag to indicate that the Q-criterion of the input vector is to be computed.
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
@ name
Definition: vtkX3D.h:225
int vtkTypeBool
Definition: vtkABI.h:69