VTK  9.5.20251125
vtkGradientFilter.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
52
53#ifndef vtkGradientFilter_h
54#define vtkGradientFilter_h
55
56#include "vtkDataSetAlgorithm.h"
57#include "vtkFiltersGeneralModule.h" // For export macro
58
59VTK_ABI_NAMESPACE_BEGIN
61
62class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
63{
64public:
66
71 void PrintSelf(ostream& os, vtkIndent indent) override;
73
76 {
77 All = 0,
78 Patch = 1,
80 };
81
85 {
86 Zero = 0,
87 NaN = 1,
90 };
91
93
99 virtual void SetInputScalars(int fieldAssociation, const char* name);
100 virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
102
104
109 vtkGetStringMacro(ResultArrayName);
110 vtkSetStringMacro(ResultArrayName);
112
114
119 vtkGetStringMacro(DivergenceArrayName);
120 vtkSetStringMacro(DivergenceArrayName);
122
124
129 vtkGetStringMacro(VorticityArrayName);
130 vtkSetStringMacro(VorticityArrayName);
132
134
139 vtkGetStringMacro(QCriterionArrayName);
140 vtkSetStringMacro(QCriterionArrayName);
142
144
157
159
166 vtkBooleanMacro(ComputeGradient, vtkTypeBool);
168
170
180
182
190 vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
192
194
205
207
211 vtkSetClampMacro(ContributingCellOption, int, 0, 2);
212 vtkGetMacro(ContributingCellOption, int);
214
216
221 vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
222 vtkGetMacro(ReplacementValueOption, int);
224
225protected:
228
231
237 virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
238 vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
239 vtkDataSet* output);
240
246 virtual int ComputeRegularGridGradient(vtkDataArray* Array, int* dims, int fieldAssociation,
247 bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output,
248 vtkUnsignedCharArray* ghosts, unsigned char hiddenGhost);
249
257
263
269
275
281
292
298
305
312
319
324 int ContributingCellOption;
325
332
333private:
334 vtkGradientFilter(const vtkGradientFilter&) = delete;
335 void operator=(const vtkGradientFilter&) = delete;
336};
337
338VTK_ABI_NAMESPACE_END
339#endif //_vtkGradientFilter_h
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
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.
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.
@ DataSetMax
Highest dimension cells in the data set.
@ Patch
Highest dimension cells in the patch of cells contributing to the 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...
@ DataTypeMax
The maximum possible value of the input array data type.
@ DataTypeMin
The minimum possible value of the input array data type.
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:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of unsigned char
int vtkTypeBool
Definition vtkABI.h:64
#define vtkDataArray
#define vtkGradientFilter