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
46#ifndef vtkGradientFilter_h
47#define vtkGradientFilter_h
48
49#include "vtkDataSetAlgorithm.h"
50#include "vtkFiltersGeneralModule.h" // For export macro
51
53
54class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
55{
56public:
58
63 void PrintSelf(ostream& os, vtkIndent indent) override;
65
68 {
69 All = 0,
70 Patch = 1,
71 DataSetMax = 2
72 };
73
77 {
78 Zero = 0,
79 NaN = 1,
80 DataTypeMin = 2,
81 DataTypeMax = 3
82 };
83
85
91 virtual void SetInputScalars(int fieldAssociation, const char* name);
92 virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
94
96
101 vtkGetStringMacro(ResultArrayName);
102 vtkSetStringMacro(ResultArrayName);
104
106
111 vtkGetStringMacro(DivergenceArrayName);
112 vtkSetStringMacro(DivergenceArrayName);
114
116
121 vtkGetStringMacro(VorticityArrayName);
122 vtkSetStringMacro(VorticityArrayName);
124
126
131 vtkGetStringMacro(QCriterionArrayName);
132 vtkSetStringMacro(QCriterionArrayName);
134
136
145 vtkGetMacro(FasterApproximation, vtkTypeBool);
146 vtkSetMacro(FasterApproximation, vtkTypeBool);
147 vtkBooleanMacro(FasterApproximation, vtkTypeBool);
149
151
156 vtkSetMacro(ComputeGradient, vtkTypeBool);
157 vtkGetMacro(ComputeGradient, vtkTypeBool);
158 vtkBooleanMacro(ComputeGradient, vtkTypeBool);
160
162
168 vtkSetMacro(ComputeDivergence, vtkTypeBool);
169 vtkGetMacro(ComputeDivergence, vtkTypeBool);
170 vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
172
174
180 vtkSetMacro(ComputeVorticity, vtkTypeBool);
181 vtkGetMacro(ComputeVorticity, vtkTypeBool);
182 vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
184
186
193 vtkSetMacro(ComputeQCriterion, vtkTypeBool);
194 vtkGetMacro(ComputeQCriterion, vtkTypeBool);
195 vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
197
199
203 vtkSetClampMacro(ContributingCellOption, int, 0, 2);
204 vtkGetMacro(ContributingCellOption, int);
206
208
213 vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
214 vtkGetMacro(ReplacementValueOption, int);
216
217protected:
220
223
229 virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
230 vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
231 vtkDataSet* output);
232
238 virtual int ComputeRegularGridGradient(vtkDataArray* Array, int* dims, int fieldAssociation,
239 bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output,
240 vtkUnsignedCharArray* ghosts, unsigned char hiddenGhost);
241
249
255
261
267
273
284
290
297
304
311
317
324
325private:
326 vtkGradientFilter(const vtkGradientFilter&) = delete;
327 void operator=(const vtkGradientFilter&) = delete;
328};
329
330#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: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.
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:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of unsigned char
int vtkTypeBool
Definition vtkABI.h:69