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
65#ifndef vtkGradientFilter_h
66#define vtkGradientFilter_h
67
68#include "vtkDataSetAlgorithm.h"
69#include "vtkFiltersGeneralModule.h" // For export macro
70
71class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
72{
73public:
75
80 void PrintSelf(ostream& os, vtkIndent indent) override;
82
85 {
86 All = 0,
87 Patch = 1,
88 DataSetMax = 2
89 };
90
94 {
95 Zero = 0,
96 NaN = 1,
97 DataTypeMin = 2,
98 DataTypeMax = 3
99 };
100
102
108 virtual void SetInputScalars(int fieldAssociation, const char* name);
109 virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
111
113
118 vtkGetStringMacro(ResultArrayName);
119 vtkSetStringMacro(ResultArrayName);
121
123
128 vtkGetStringMacro(DivergenceArrayName);
129 vtkSetStringMacro(DivergenceArrayName);
131
133
138 vtkGetStringMacro(VorticityArrayName);
139 vtkSetStringMacro(VorticityArrayName);
141
143
148 vtkGetStringMacro(QCriterionArrayName);
149 vtkSetStringMacro(QCriterionArrayName);
151
153
162 vtkGetMacro(FasterApproximation, vtkTypeBool);
163 vtkSetMacro(FasterApproximation, vtkTypeBool);
164 vtkBooleanMacro(FasterApproximation, vtkTypeBool);
166
168
173 vtkSetMacro(ComputeGradient, vtkTypeBool);
174 vtkGetMacro(ComputeGradient, vtkTypeBool);
175 vtkBooleanMacro(ComputeGradient, vtkTypeBool);
177
179
185 vtkSetMacro(ComputeDivergence, vtkTypeBool);
186 vtkGetMacro(ComputeDivergence, vtkTypeBool);
187 vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
189
191
197 vtkSetMacro(ComputeVorticity, vtkTypeBool);
198 vtkGetMacro(ComputeVorticity, vtkTypeBool);
199 vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
201
203
210 vtkSetMacro(ComputeQCriterion, vtkTypeBool);
211 vtkGetMacro(ComputeQCriterion, vtkTypeBool);
212 vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
214
216
220 vtkSetClampMacro(ContributingCellOption, int, 0, 2);
221 vtkGetMacro(ContributingCellOption, int);
223
225
230 vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
231 vtkGetMacro(ReplacementValueOption, int);
233
234protected:
237
240
246 virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
247 vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
248 vtkDataSet* output);
249
255 virtual int ComputeRegularGridGradient(vtkDataArray* Array, int fieldAssociation,
256 bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output);
257
265
271
277
283
289
300
306
313
320
327
333
340
341private:
342 vtkGradientFilter(const vtkGradientFilter&) = delete;
343 void operator=(const vtkGradientFilter&) = delete;
344};
345
346#endif //_vtkGradientFilter_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
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.
static vtkGradientFilter * New()
Standard methods for instantiation, obtaining type information, and printing.
ContributingCellEnum
Options to choose what cells contribute to the gradient calculation.
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:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
int vtkTypeBool
Definition: vtkABI.h:69