VTK  9.0.3
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 
44 #ifndef vtkGradientFilter_h
45 #define vtkGradientFilter_h
46 
47 #include "vtkDataSetAlgorithm.h"
48 #include "vtkFiltersGeneralModule.h" // For export macro
49 
50 class VTKFILTERSGENERAL_EXPORT vtkGradientFilter : public vtkDataSetAlgorithm
51 {
52 public:
54  void PrintSelf(ostream& os, vtkIndent indent) override;
55 
58  {
59  All = 0,
60  Patch = 1,
61  DataSetMax = 2
62  };
63 
67  {
68  Zero = 0,
69  NaN = 1,
70  DataTypeMin = 2,
71  DataTypeMax = 3
72  };
73 
75 
77 
83  virtual void SetInputScalars(int fieldAssociation, const char* name);
84  virtual void SetInputScalars(int fieldAssociation, int fieldAttributeType);
86 
88 
93  vtkGetStringMacro(ResultArrayName);
94  vtkSetStringMacro(ResultArrayName);
96 
98 
103  vtkGetStringMacro(DivergenceArrayName);
104  vtkSetStringMacro(DivergenceArrayName);
106 
108 
113  vtkGetStringMacro(VorticityArrayName);
114  vtkSetStringMacro(VorticityArrayName);
116 
118 
123  vtkGetStringMacro(QCriterionArrayName);
124  vtkSetStringMacro(QCriterionArrayName);
126 
128 
137  vtkGetMacro(FasterApproximation, vtkTypeBool);
138  vtkSetMacro(FasterApproximation, vtkTypeBool);
139  vtkBooleanMacro(FasterApproximation, vtkTypeBool);
141 
143 
148  vtkSetMacro(ComputeGradient, vtkTypeBool);
149  vtkGetMacro(ComputeGradient, vtkTypeBool);
150  vtkBooleanMacro(ComputeGradient, vtkTypeBool);
152 
154 
160  vtkSetMacro(ComputeDivergence, vtkTypeBool);
161  vtkGetMacro(ComputeDivergence, vtkTypeBool);
162  vtkBooleanMacro(ComputeDivergence, vtkTypeBool);
164 
166 
172  vtkSetMacro(ComputeVorticity, vtkTypeBool);
173  vtkGetMacro(ComputeVorticity, vtkTypeBool);
174  vtkBooleanMacro(ComputeVorticity, vtkTypeBool);
176 
178 
185  vtkSetMacro(ComputeQCriterion, vtkTypeBool);
186  vtkGetMacro(ComputeQCriterion, vtkTypeBool);
187  vtkBooleanMacro(ComputeQCriterion, vtkTypeBool);
189 
191 
195  vtkSetClampMacro(ContributingCellOption, int, 0, 2);
196  vtkGetMacro(ContributingCellOption, int);
198 
200 
205  vtkSetClampMacro(ReplacementValueOption, int, 0, 3);
206  vtkGetMacro(ReplacementValueOption, int);
208 
209 protected:
211  ~vtkGradientFilter() override;
212 
215 
221  virtual int ComputeUnstructuredGridGradient(vtkDataArray* Array, int fieldAssociation,
222  vtkDataSet* input, bool computeVorticity, bool computeQCriterion, bool computeDivergence,
223  vtkDataSet* output);
224 
230  virtual int ComputeRegularGridGradient(vtkDataArray* Array, int fieldAssociation,
231  bool computeVorticity, bool computeQCriterion, bool computeDivergence, vtkDataSet* output);
232 
240 
246 
252 
258 
264 
275 
281 
288 
295 
302 
308 
315 
316 private:
317  vtkGradientFilter(const vtkGradientFilter&) = delete;
318  void operator=(const vtkGradientFilter&) = delete;
319 };
320 
321 #endif //_vtkGradientFilter_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
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)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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()
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:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
@ name
Definition: vtkX3D.h:225
int vtkTypeBool
Definition: vtkABI.h:69