VTK  9.2.6
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
49#ifndef vtkGradientFilter_h
50#define vtkGradientFilter_h
51
52#include "vtkDataSetAlgorithm.h"
53#include "vtkFiltersGeneralModule.h" // For export macro
54
56
57class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
58{
59public:
61
66 void PrintSelf(ostream& os, vtkIndent indent) override;
68
71 {
72 All = 0,
73 Patch = 1,
74 DataSetMax = 2
75 };
76
80 {
81 Zero = 0,
82 NaN = 1,
83 DataTypeMin = 2,
84 DataTypeMax = 3
85 };
86
88
94 virtual void SetInputScalars(int fieldAssociation, const char* name);
95 virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
97
99
104 vtkGetStringMacro(ResultArrayName);
105 vtkSetStringMacro(ResultArrayName);
107
109
114 vtkGetStringMacro(DivergenceArrayName);
115 vtkSetStringMacro(DivergenceArrayName);
117
119
124 vtkGetStringMacro(VorticityArrayName);
125 vtkSetStringMacro(VorticityArrayName);
127
129
134 vtkGetStringMacro(QCriterionArrayName);
135 vtkSetStringMacro(QCriterionArrayName);
137
139
148 vtkGetMacro(FasterApproximation, vtkTypeBool);
149 vtkSetMacro(FasterApproximation, vtkTypeBool);
150 vtkBooleanMacro(FasterApproximation, vtkTypeBool);
152
154
159 vtkSetMacro(ComputeGradient, vtkTypeBool);
160 vtkGetMacro(ComputeGradient, vtkTypeBool);
161 vtkBooleanMacro(ComputeGradient, vtkTypeBool);
163
165
171 vtkSetMacro(ComputeDivergence, vtkTypeBool);
172 vtkGetMacro(ComputeDivergence, vtkTypeBool);
173 vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
175
177
183 vtkSetMacro(ComputeVorticity, vtkTypeBool);
184 vtkGetMacro(ComputeVorticity, vtkTypeBool);
185 vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
187
189
196 vtkSetMacro(ComputeQCriterion, vtkTypeBool);
197 vtkGetMacro(ComputeQCriterion, vtkTypeBool);
198 vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
200
202
206 vtkSetClampMacro(ContributingCellOption, int, 0, 2);
207 vtkGetMacro(ContributingCellOption, int);
209
211
216 vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
217 vtkGetMacro(ReplacementValueOption, int);
219
220protected:
223
226
232 virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
233 vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
234 vtkDataSet* output);
235
241 virtual int ComputeRegularGridGradient(vtkDataArray* Array, int* dims, int fieldAssociation,
242 bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output,
243 vtkUnsignedCharArray* ghosts, unsigned char hiddenGhost);
244
252
258
264
270
276
287
293
300
307
314
320
327
328private:
329 vtkGradientFilter(const vtkGradientFilter&) = delete;
330 void operator=(const vtkGradientFilter&) = delete;
331};
332
333#endif //_vtkGradientFilter_h
abstract superclass for arrays of numeric data
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition vtkDataSet.h:63
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 int ComputeRegularGridGradient(vtkDataArray *Array, int *dims, int fieldAssociation, bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet *output, vtkUnsignedCharArray *ghosts, unsigned char hiddenGhost)
Compute the gradients for either a vtkImageData, vtkRectilinearGrid or a vtkStructuredGrid.
virtual void SetInputScalars(int fieldAssociation, const char *name)
These are basically a convenience method that calls SetInputArrayToProcess to set the array used as t...
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:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of unsigned char
int vtkTypeBool
Definition vtkABI.h:69