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 
68 #ifndef vtkGradientFilter_h
69 #define vtkGradientFilter_h
70 
71 #include "vtkDataSetAlgorithm.h"
72 #include "vtkFiltersGeneralModule.h" // For export macro
73 
75 
76 class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
77 {
78 public:
80 
85  void PrintSelf(ostream& os, vtkIndent indent) override;
87 
90  {
91  All = 0,
92  Patch = 1,
93  DataSetMax = 2
94  };
95 
99  {
100  Zero = 0,
101  NaN = 1,
102  DataTypeMin = 2,
103  DataTypeMax = 3
104  };
105 
107 
113  virtual void SetInputScalars(int fieldAssociation, const char* name);
114  virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
116 
118 
123  vtkGetStringMacro(ResultArrayName);
124  vtkSetStringMacro(ResultArrayName);
126 
128 
133  vtkGetStringMacro(DivergenceArrayName);
134  vtkSetStringMacro(DivergenceArrayName);
136 
138 
143  vtkGetStringMacro(VorticityArrayName);
144  vtkSetStringMacro(VorticityArrayName);
146 
148 
153  vtkGetStringMacro(QCriterionArrayName);
154  vtkSetStringMacro(QCriterionArrayName);
156 
158 
167  vtkGetMacro(FasterApproximation, vtkTypeBool);
168  vtkSetMacro(FasterApproximation, vtkTypeBool);
169  vtkBooleanMacro(FasterApproximation, vtkTypeBool);
171 
173 
178  vtkSetMacro(ComputeGradient, vtkTypeBool);
179  vtkGetMacro(ComputeGradient, vtkTypeBool);
180  vtkBooleanMacro(ComputeGradient, vtkTypeBool);
182 
184 
190  vtkSetMacro(ComputeDivergence, vtkTypeBool);
191  vtkGetMacro(ComputeDivergence, vtkTypeBool);
192  vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
194 
196 
202  vtkSetMacro(ComputeVorticity, vtkTypeBool);
203  vtkGetMacro(ComputeVorticity, vtkTypeBool);
204  vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
206 
208 
215  vtkSetMacro(ComputeQCriterion, vtkTypeBool);
216  vtkGetMacro(ComputeQCriterion, vtkTypeBool);
217  vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
219 
221 
225  vtkSetClampMacro(ContributingCellOption, int, 0, 2);
226  vtkGetMacro(ContributingCellOption, int);
228 
230 
235  vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
236  vtkGetMacro(ReplacementValueOption, int);
238 
239 protected:
241  ~vtkGradientFilter() override;
242 
245 
251  virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
252  vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
253  vtkDataSet* output);
254 
260  virtual int ComputeRegularGridGradient(vtkDataArray* Array, int* dims, int fieldAssociation,
261  bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output,
262  vtkUnsignedCharArray* ghosts, unsigned char hiddenGhost);
263 
271 
277 
283 
289 
295 
306 
312 
319 
326 
333 
339 
346 
347 private:
348  vtkGradientFilter(const vtkGradientFilter&) = delete;
349  void operator=(const vtkGradientFilter&) = delete;
350 };
351 
352 #endif //_vtkGradientFilter_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:172
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.
ContributingCellEnum
Options to choose what cells contribute to the gradient calculation.
static vtkGradientFilter * New()
Standard methods for instantiation, obtaining type information, and printing.
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:119
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of unsigned char
@ name
Definition: vtkX3D.h:225
int vtkTypeBool
Definition: vtkABI.h:69