VTK
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkGradientFilter Class Reference

A general filter for gradient estimation. More...

#include <vtkGradientFilter.h>

Inherits vtkDataSetAlgorithm.

Public Types

typedef vtkDataSetAlgorithm Superclass
 

Public Member Functions

virtual int IsA (const char *type)
 
vtkGradientFilterNewInstance () const
 
virtual void PrintSelf (ostream &os, vtkIndent indent)
 
virtual void SetInputScalars (int fieldAssociation, const char *name)
 
virtual void SetInputScalars (int fieldAssociation, int fieldAttributeType)
 
virtual char * GetResultArrayName ()
 
virtual void SetResultArrayName (const char *)
 
virtual char * GetVorticityArrayName ()
 
virtual void SetVorticityArrayName (const char *)
 
virtual char * GetQCriterionArrayName ()
 
virtual void SetQCriterionArrayName (const char *)
 
virtual int GetFasterApproximation ()
 
virtual void SetFasterApproximation (int)
 
virtual void FasterApproximationOn ()
 
virtual void FasterApproximationOff ()
 
virtual void SetComputeVorticity (int)
 
virtual int GetComputeVorticity ()
 
virtual void ComputeVorticityOn ()
 
virtual void ComputeVorticityOff ()
 
virtual void SetComputeQCriterion (int)
 
virtual int GetComputeQCriterion ()
 
virtual void ComputeQCriterionOn ()
 
virtual void ComputeQCriterionOff ()
 

Static Public Member Functions

static int IsTypeOf (const char *type)
 
static vtkGradientFilterSafeDownCast (vtkObjectBase *o)
 
static vtkGradientFilterNew ()
 

Protected Member Functions

virtual vtkObjectBase * NewInstanceInternal () const
 
 vtkGradientFilter ()
 
 ~vtkGradientFilter ()
 
virtual int RequestUpdateExtent (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 
virtual int RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 
virtual int ComputeUnstructuredGridGradient (vtkDataArray *Array, int fieldAssociation, vtkDataSet *input, bool computeVorticity, bool computeQCriterion, vtkDataSet *output)
 
virtual int ComputeRegularGridGradient (vtkDataArray *Array, int fieldAssociation, bool computeVorticity, bool computeQCriterion, vtkDataSet *output)
 

Protected Attributes

char * ResultArrayName
 
char * VorticityArrayName
 
char * QCriterionArrayName
 
int FasterApproximation
 
int ComputeQCriterion
 
int ComputeVorticity
 

Detailed Description

A general filter for gradient estimation.

Estimates the gradient of a field in a data set. The gradient calculation is dependent on the input dataset type. The created gradient array is of the same type as the array it is calculated from (e.g. point data or cell data) as well as data type (e.g. float, double). At the boundary the gradient is not central differencing. The output array has 3*number of components of the input data array. The ordering for the output tuple will be {du/dx, du/dy, du/dz, dv/dx, dv/dy, dv/dz, dw/dx, dw/dy, dw/dz} for an input array {u, v, w}. There are also the options to additionally compute the vorticity and Q criterion of a vector field.

Tests:
vtkGradientFilter (Tests)

Definition at line 44 of file vtkGradientFilter.h.

Member Typedef Documentation

typedef vtkDataSetAlgorithm vtkGradientFilter::Superclass

Definition at line 47 of file vtkGradientFilter.h.

Constructor & Destructor Documentation

vtkGradientFilter::vtkGradientFilter ( )
protected
vtkGradientFilter::~vtkGradientFilter ( )
protected

Member Function Documentation

static int vtkGradientFilter::IsTypeOf ( const char *  type)
static
virtual int vtkGradientFilter::IsA ( const char *  type)
virtual
static vtkGradientFilter* vtkGradientFilter::SafeDownCast ( vtkObjectBase *  o)
static
virtual vtkObjectBase* vtkGradientFilter::NewInstanceInternal ( ) const
protectedvirtual
vtkGradientFilter* vtkGradientFilter::NewInstance ( ) const
virtual void vtkGradientFilter::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
virtual
static vtkGradientFilter* vtkGradientFilter::New ( )
static
virtual void vtkGradientFilter::SetInputScalars ( int  fieldAssociation,
const char *  name 
)
virtual

These are basically a convenience method that calls SetInputArrayToProcess to set the array used as the input scalars. The fieldAssociation comes from the vtkDataObject::FieldAssocations enum. The fieldAttributeType comes from the vtkDataSetAttributes::AttributeTypes enum.

virtual void vtkGradientFilter::SetInputScalars ( int  fieldAssociation,
int  fieldAttributeType 
)
virtual

These are basically a convenience method that calls SetInputArrayToProcess to set the array used as the input scalars. The fieldAssociation comes from the vtkDataObject::FieldAssocations enum. The fieldAttributeType comes from the vtkDataSetAttributes::AttributeTypes enum.

virtual char* vtkGradientFilter::GetResultArrayName ( )
virtual

Get/Set the name of the gradient array to create. If NULL (the default) then the output array will be named "Gradients".

virtual void vtkGradientFilter::SetResultArrayName ( const char *  )
virtual

Get/Set the name of the gradient array to create. If NULL (the default) then the output array will be named "Gradients".

virtual char* vtkGradientFilter::GetVorticityArrayName ( )
virtual

Get/Set the name of the vorticity array to create. This is only used if ComputeVorticity is non-zero. If NULL (the default) then the output array will be named "Vorticity".

virtual void vtkGradientFilter::SetVorticityArrayName ( const char *  )
virtual

Get/Set the name of the vorticity array to create. This is only used if ComputeVorticity is non-zero. If NULL (the default) then the output array will be named "Vorticity".

virtual char* vtkGradientFilter::GetQCriterionArrayName ( )
virtual

Get/Set the name of the Q criterion array to create. This is only used if ComputeQCriterion is non-zero. If NULL (the default) then the output array will be named "Q-criterion".

virtual void vtkGradientFilter::SetQCriterionArrayName ( const char *  )
virtual

Get/Set the name of the Q criterion array to create. This is only used if ComputeQCriterion is non-zero. If NULL (the default) then the output array will be named "Q-criterion".

virtual int vtkGradientFilter::GetFasterApproximation ( )
virtual

When this flag is on (default is off), the gradient filter will provide a less accurate (but close) algorithm that performs fewer derivative calculations (and is therefore faster). The error contains some smoothing of the output data and some possible errors on the boundary. This parameter has no effect when performing the gradient of cell data. This only applies if the input grid is a vtkUnstructuredGrid or a vtkPolyData.

virtual void vtkGradientFilter::SetFasterApproximation ( int  )
virtual

When this flag is on (default is off), the gradient filter will provide a less accurate (but close) algorithm that performs fewer derivative calculations (and is therefore faster). The error contains some smoothing of the output data and some possible errors on the boundary. This parameter has no effect when performing the gradient of cell data. This only applies if the input grid is a vtkUnstructuredGrid or a vtkPolyData.

virtual void vtkGradientFilter::FasterApproximationOn ( )
virtual

When this flag is on (default is off), the gradient filter will provide a less accurate (but close) algorithm that performs fewer derivative calculations (and is therefore faster). The error contains some smoothing of the output data and some possible errors on the boundary. This parameter has no effect when performing the gradient of cell data. This only applies if the input grid is a vtkUnstructuredGrid or a vtkPolyData.

virtual void vtkGradientFilter::FasterApproximationOff ( )
virtual

When this flag is on (default is off), the gradient filter will provide a less accurate (but close) algorithm that performs fewer derivative calculations (and is therefore faster). The error contains some smoothing of the output data and some possible errors on the boundary. This parameter has no effect when performing the gradient of cell data. This only applies if the input grid is a vtkUnstructuredGrid or a vtkPolyData.

virtual void vtkGradientFilter::SetComputeVorticity ( int  )
virtual

Set the resultant array to be vorticity/curl of the input array. The input array must have 3 components.

virtual int vtkGradientFilter::GetComputeVorticity ( )
virtual

Set the resultant array to be vorticity/curl of the input array. The input array must have 3 components.

virtual void vtkGradientFilter::ComputeVorticityOn ( )
virtual

Set the resultant array to be vorticity/curl of the input array. The input array must have 3 components.

virtual void vtkGradientFilter::ComputeVorticityOff ( )
virtual

Set the resultant array to be vorticity/curl of the input array. The input array must have 3 components.

virtual void vtkGradientFilter::SetComputeQCriterion ( int  )
virtual

Add Q-criterion to the output field data. The name of the array will be "Q-Criterion" and will be the same type as the input array. The input array must have 3 components in order to compute this. Note that Q-citerion is a balance of the rate of vorticity and the rate of strain.

virtual int vtkGradientFilter::GetComputeQCriterion ( )
virtual

Add Q-criterion to the output field data. The name of the array will be "Q-Criterion" and will be the same type as the input array. The input array must have 3 components in order to compute this. Note that Q-citerion is a balance of the rate of vorticity and the rate of strain.

virtual void vtkGradientFilter::ComputeQCriterionOn ( )
virtual

Add Q-criterion to the output field data. The name of the array will be "Q-Criterion" and will be the same type as the input array. The input array must have 3 components in order to compute this. Note that Q-citerion is a balance of the rate of vorticity and the rate of strain.

virtual void vtkGradientFilter::ComputeQCriterionOff ( )
virtual

Add Q-criterion to the output field data. The name of the array will be "Q-Criterion" and will be the same type as the input array. The input array must have 3 components in order to compute this. Note that Q-citerion is a balance of the rate of vorticity and the rate of strain.

virtual int vtkGradientFilter::RequestUpdateExtent ( vtkInformation *  ,
vtkInformationVector **  ,
vtkInformationVector *   
)
protectedvirtual
virtual int vtkGradientFilter::RequestData ( vtkInformation *  ,
vtkInformationVector **  ,
vtkInformationVector *   
)
protectedvirtual
virtual int vtkGradientFilter::ComputeUnstructuredGridGradient ( vtkDataArray *  Array,
int  fieldAssociation,
vtkDataSet *  input,
bool  computeVorticity,
bool  computeQCriterion,
vtkDataSet *  output 
)
protectedvirtual

Compute the gradients for grids that are not a vtkImageData, vtkRectilinearGrid, or vtkStructuredGrid. Returns non-zero if the operation was successful.

virtual int vtkGradientFilter::ComputeRegularGridGradient ( vtkDataArray *  Array,
int  fieldAssociation,
bool  computeVorticity,
bool  computeQCriterion,
vtkDataSet *  output 
)
protectedvirtual

Compute the gradients for either a vtkImageData, vtkRectilinearGrid or a vtkStructuredGrid. Computes the gradient using finite differences. Returns non-zero if the operation was successful.

Member Data Documentation

char* vtkGradientFilter::ResultArrayName
protected

If non-null then it contains the name of the outputted gradient array. By derault it is "Gradients".

Definition at line 147 of file vtkGradientFilter.h.

char* vtkGradientFilter::VorticityArrayName
protected

If non-null then it contains the name of the outputted vorticity array. By derault it is "Vorticity".

Definition at line 151 of file vtkGradientFilter.h.

char* vtkGradientFilter::QCriterionArrayName
protected

If non-null then it contains the name of the outputted Q criterion array. By derault it is "Q-criterion".

Definition at line 155 of file vtkGradientFilter.h.

int vtkGradientFilter::FasterApproximation
protected

When this flag is on (default is off), the gradient filter will provide a less accurate (but close) algorithm that performs fewer derivative calculations (and is therefore faster). The error contains some smoothing of the output data and some possible errors on the boundary. This parameter has no effect when performing the gradient of cell data. This only applies if the input grid is a vtkUnstructuredGrid or a vtkPolyData.

Definition at line 164 of file vtkGradientFilter.h.

int vtkGradientFilter::ComputeQCriterion
protected

Flag to indicate that the Q-criterion of the input vector is to be computed. The input array to be processed must have 3 components. By default ComputeQCriterion is off.

Definition at line 169 of file vtkGradientFilter.h.

int vtkGradientFilter::ComputeVorticity
protected

Flag to indicate that vorticity/curl of the input vector is to be computed. The input array to be processed must have 3 components. By default ComputeVorticity is off.

Definition at line 174 of file vtkGradientFilter.h.


The documentation for this class was generated from the following file: