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
55#ifndef vtkGradientFilter_h
56#define vtkGradientFilter_h
57
58#include "vtkDataSetAlgorithm.h"
59#include "vtkFiltersGeneralModule.h" // For export macro
60
61class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
62{
63public:
65
70 void PrintSelf(ostream& os, vtkIndent indent) override;
72
75 {
76 All = 0,
77 Patch = 1,
78 DataSetMax = 2
79 };
80
84 {
85 Zero = 0,
86 NaN = 1,
87 DataTypeMin = 2,
88 DataTypeMax = 3
89 };
90
92
98 virtual void SetInputScalars(int fieldAssociation, const char* name);
99 virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
101
103
108 vtkGetStringMacro(ResultArrayName);
109 vtkSetStringMacro(ResultArrayName);
111
113
118 vtkGetStringMacro(DivergenceArrayName);
119 vtkSetStringMacro(DivergenceArrayName);
121
123
128 vtkGetStringMacro(VorticityArrayName);
129 vtkSetStringMacro(VorticityArrayName);
131
133
138 vtkGetStringMacro(QCriterionArrayName);
139 vtkSetStringMacro(QCriterionArrayName);
141
143
152 vtkGetMacro(FasterApproximation, vtkTypeBool);
153 vtkSetMacro(FasterApproximation, vtkTypeBool);
154 vtkBooleanMacro(FasterApproximation, vtkTypeBool);
156
158
163 vtkSetMacro(ComputeGradient, vtkTypeBool);
164 vtkGetMacro(ComputeGradient, vtkTypeBool);
165 vtkBooleanMacro(ComputeGradient, vtkTypeBool);
167
169
175 vtkSetMacro(ComputeDivergence, vtkTypeBool);
176 vtkGetMacro(ComputeDivergence, vtkTypeBool);
177 vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
179
181
187 vtkSetMacro(ComputeVorticity, vtkTypeBool);
188 vtkGetMacro(ComputeVorticity, vtkTypeBool);
189 vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
191
193
200 vtkSetMacro(ComputeQCriterion, vtkTypeBool);
201 vtkGetMacro(ComputeQCriterion, vtkTypeBool);
202 vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
204
206
210 vtkSetClampMacro(ContributingCellOption, int, 0, 2);
211 vtkGetMacro(ContributingCellOption, int);
213
215
220 vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
221 vtkGetMacro(ReplacementValueOption, int);
223
224protected:
227
230
236 virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
237 vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
238 vtkDataSet* output);
239
245 virtual int ComputeRegularGridGradient(vtkDataArray* Array, int fieldAssociation,
246 bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output);
247
255
261
267
273
279
290
296
303
310
317
323
330
331private:
332 vtkGradientFilter(const vtkGradientFilter&) = delete;
333 void operator=(const vtkGradientFilter&) = delete;
334};
335
336#endif //_vtkGradientFilter_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:66
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:43
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
@ name
Definition: vtkX3D.h:225
int vtkTypeBool
Definition: vtkABI.h:69