VTK-m  2.2
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
vtkm::filter::Filter Class Referenceabstract

Base class for all filters. More...

#include <Filter.h>

Inheritance diagram for vtkm::filter::Filter:
vtkm::filter::clean_grid::CleanGrid vtkm::filter::connected_components::CellSetConnectivity vtkm::filter::connected_components::ImageConnectivity vtkm::filter::contour::AbstractContour vtkm::filter::contour::ClipWithField vtkm::filter::contour::ClipWithImplicitFunction vtkm::filter::contour::MIRFilter vtkm::filter::density_estimate::ContinuousScatterPlot vtkm::filter::density_estimate::Entropy vtkm::filter::density_estimate::Histogram vtkm::filter::density_estimate::NDEntropy vtkm::filter::density_estimate::NDHistogram vtkm::filter::density_estimate::ParticleDensityBase vtkm::filter::density_estimate::Statistics vtkm::filter::entity_extraction::ExternalFaces vtkm::filter::entity_extraction::ExtractGeometry vtkm::filter::entity_extraction::ExtractPoints vtkm::filter::entity_extraction::ExtractStructured vtkm::filter::entity_extraction::GhostCellRemove vtkm::filter::entity_extraction::Mask vtkm::filter::entity_extraction::MaskPoints vtkm::filter::entity_extraction::Threshold vtkm::filter::entity_extraction::ThresholdPoints vtkm::filter::field_conversion::CellAverage vtkm::filter::field_conversion::PointAverage vtkm::filter::field_transform::CompositeVectors vtkm::filter::field_transform::CylindricalCoordinateTransform vtkm::filter::field_transform::FieldToColors vtkm::filter::field_transform::GenerateIds vtkm::filter::field_transform::LogValues vtkm::filter::field_transform::PointElevation vtkm::filter::field_transform::PointTransform vtkm::filter::field_transform::SphericalCoordinateTransform vtkm::filter::field_transform::Warp vtkm::filter::FilterField vtkm::filter::flow::FilterParticleAdvection vtkm::filter::flow::Lagrangian vtkm::filter::flow::LagrangianStructures vtkm::filter::flow::StreamSurface vtkm::filter::geometry_refinement::ConvertToPointCloud vtkm::filter::geometry_refinement::Shrink vtkm::filter::geometry_refinement::SplitSharpEdges vtkm::filter::geometry_refinement::Tetrahedralize vtkm::filter::geometry_refinement::Triangulate vtkm::filter::geometry_refinement::Tube vtkm::filter::geometry_refinement::VertexClustering vtkm::filter::image_processing::ComputeMoments vtkm::filter::image_processing::ImageDifference vtkm::filter::image_processing::ImageMedian vtkm::filter::mesh_info::CellMeasures vtkm::filter::mesh_info::GhostCellClassify vtkm::filter::mesh_info::MeshQuality vtkm::filter::mesh_info::MeshQualityArea vtkm::filter::mesh_info::MeshQualityAspectGamma vtkm::filter::mesh_info::MeshQualityAspectRatio vtkm::filter::mesh_info::MeshQualityCondition vtkm::filter::mesh_info::MeshQualityDiagonalRatio vtkm::filter::mesh_info::MeshQualityDimension vtkm::filter::mesh_info::MeshQualityJacobian vtkm::filter::mesh_info::MeshQualityMaxAngle vtkm::filter::mesh_info::MeshQualityMaxDiagonal vtkm::filter::mesh_info::MeshQualityMinAngle vtkm::filter::mesh_info::MeshQualityMinDiagonal vtkm::filter::mesh_info::MeshQualityOddy vtkm::filter::mesh_info::MeshQualityRelativeSizeSquared vtkm::filter::mesh_info::MeshQualityScaledJacobian vtkm::filter::mesh_info::MeshQualityShape vtkm::filter::mesh_info::MeshQualityShapeAndSize vtkm::filter::mesh_info::MeshQualityShear vtkm::filter::mesh_info::MeshQualitySkew vtkm::filter::mesh_info::MeshQualityStretch vtkm::filter::mesh_info::MeshQualityTaper vtkm::filter::mesh_info::MeshQualityVolume vtkm::filter::mesh_info::MeshQualityWarpage vtkm::filter::multi_block::AmrArrays vtkm::filter::multi_block::MergeDataSets vtkm::filter::resampling::HistSampling vtkm::filter::resampling::Probe vtkm::filter::scalar_topology::ContourTreeAugmented vtkm::filter::scalar_topology::ContourTreeMesh2D vtkm::filter::scalar_topology::ContourTreeMesh3D vtkm::filter::scalar_topology::ContourTreeUniformDistributed vtkm::filter::scalar_topology::DistributedBranchDecompositionFilter vtkm::filter::scalar_topology::SelectTopVolumeContoursFilter vtkm::filter::uncertainty::ContourUncertainUniform vtkm::filter::uncertainty::ContourUncertainUniformMonteCarlo vtkm::filter::vector_analysis::CrossProduct vtkm::filter::vector_analysis::DotProduct vtkm::filter::vector_analysis::Gradient vtkm::filter::vector_analysis::SurfaceNormals vtkm::filter::vector_analysis::VectorMagnitude vtkm::filter::zfp::ZFPCompressor1D vtkm::filter::zfp::ZFPCompressor2D vtkm::filter::zfp::ZFPCompressor3D vtkm::filter::zfp::ZFPDecompressor1D vtkm::filter::zfp::ZFPDecompressor2D vtkm::filter::zfp::ZFPDecompressor3D

Classes

struct  ScalarToVec
 

Public Member Functions

 Filter ()
 
virtual ~Filter ()
 
vtkm::cont::DataSet Execute (const vtkm::cont::DataSet &input)
 Executes the filter on the input and produces a result dataset. More...
 
vtkm::cont::PartitionedDataSet Execute (const vtkm::cont::PartitionedDataSet &input)
 Executes the filter on the input PartitionedDataSet and produces a result PartitionedDataSet. More...
 
void SetFieldsToPass (const vtkm::filter::FieldSelection &fieldsToPass)
 Specify which fields get passed from input to output. More...
 
void SetFieldsToPass (vtkm::filter::FieldSelection &&fieldsToPass)
 Specify which fields get passed from input to output. More...
 
void SetFieldsToPass (const vtkm::filter::FieldSelection &fieldsToPass, vtkm::filter::FieldSelection::Mode mode)
 
void SetFieldsToPass (std::initializer_list< std::string > fields, vtkm::filter::FieldSelection::Mode mode=vtkm::filter::FieldSelection::Mode::Select)
 Specify which fields get passed from input to output. More...
 
void SetFieldsToPass (std::initializer_list< std::pair< std::string, vtkm::cont::Field::Association >> fields, vtkm::filter::FieldSelection::Mode mode=vtkm::filter::FieldSelection::Mode::Select)
 Specify which fields get passed from input to output. More...
 
void SetFieldsToPass (const std::string &fieldname, vtkm::cont::Field::Association association, vtkm::filter::FieldSelection::Mode mode=vtkm::filter::FieldSelection::Mode::Select)
 Specify which fields get passed from input to output. More...
 
void SetFieldsToPass (const std::string &fieldname, vtkm::filter::FieldSelection::Mode mode)
 Specify which fields get passed from input to output. More...
 
const vtkm::filter::FieldSelectionGetFieldsToPass () const
 Specify which fields get passed from input to output. More...
 
vtkm::filter::FieldSelectionGetFieldsToPass ()
 Specify which fields get passed from input to output. More...
 
void SetPassCoordinateSystems (bool flag)
 Specify whether to always pass coordinate systems. More...
 
bool GetPassCoordinateSystems () const
 Specify whether to always pass coordinate systems. More...
 
void SetOutputFieldName (const std::string &name)
 Specifies the name of the output field generated. More...
 
const std::string & GetOutputFieldName () const
 Specifies the name of the output field generated. More...
 
void SetActiveField (const std::string &name, vtkm::cont::Field::Association association=vtkm::cont::Field::Association::Any)
 Specifies a field to operate on. More...
 
void SetActiveField (vtkm::IdComponent index, const std::string &name, vtkm::cont::Field::Association association=vtkm::cont::Field::Association::Any)
 Specifies a field to operate on. More...
 
const std::string & GetActiveFieldName (vtkm::IdComponent index=0) const
 Specifies a field to operate on. More...
 
vtkm::cont::Field::Association GetActiveFieldAssociation (vtkm::IdComponent index=0) const
 Specifies a field to operate on. More...
 
void SetActiveCoordinateSystem (vtkm::Id coord_idx)
 Specifies the coordinate system index to make active to use when processing the input vtkm::cont::DataSet. More...
 
void SetActiveCoordinateSystem (vtkm::IdComponent index, vtkm::Id coord_idx)
 Specifies the coordinate system index to make active to use when processing the input vtkm::cont::DataSet. More...
 
vtkm::Id GetActiveCoordinateSystemIndex (vtkm::IdComponent index=0) const
 Specifies the coordinate system index to make active to use when processing the input vtkm::cont::DataSet. More...
 
void SetUseCoordinateSystemAsField (bool val)
 Specifies whether to use point coordinates as the input field. More...
 
void SetUseCoordinateSystemAsField (vtkm::IdComponent index, bool val)
 Specifies whether to use point coordinates as the input field. More...
 
bool GetUseCoordinateSystemAsField (vtkm::IdComponent index=0) const
 Specifies whether to use point coordinates as the input field. More...
 
vtkm::IdComponent GetNumberOfActiveFields () const
 Return the number of active fields currently set. More...
 
virtual bool CanThread () const
 Returns whether the filter can execute on partitions in concurrent threads. More...
 
void SetThreadsPerCPU (vtkm::Id numThreads)
 
void SetThreadsPerGPU (vtkm::Id numThreads)
 
vtkm::Id GetThreadsPerCPU () const
 
vtkm::Id GetThreadsPerGPU () const
 
bool GetRunMultiThreadedFilter () const
 
void SetRunMultiThreadedFilter (bool val)
 
void SetInvoker (vtkm::cont::Invoker inv)
 Specify the vtkm::cont::Invoker to be used to execute worklets by this filter instance. More...
 

Protected Member Functions

vtkm::cont::DataSet CreateResult (const vtkm::cont::DataSet &inDataSet) const
 Create the output data set for DoExecute. More...
 
vtkm::cont::DataSet CreateResultField (const vtkm::cont::DataSet &inDataSet, const vtkm::cont::Field &resultField) const
 Create the output data set for DoExecute More...
 
vtkm::cont::DataSet CreateResultField (const vtkm::cont::DataSet &inDataSet, const std::string &resultFieldName, vtkm::cont::Field::Association resultFieldAssociation, const vtkm::cont::UnknownArrayHandle &resultFieldArray) const
 Create the output data set for DoExecute More...
 
vtkm::cont::DataSet CreateResultFieldPoint (const vtkm::cont::DataSet &inDataSet, const std::string &resultFieldName, const vtkm::cont::UnknownArrayHandle &resultFieldArray) const
 Create the output data set for DoExecute More...
 
vtkm::cont::DataSet CreateResultFieldCell (const vtkm::cont::DataSet &inDataSet, const std::string &resultFieldName, const vtkm::cont::UnknownArrayHandle &resultFieldArray) const
 Create the output data set for DoExecute More...
 
vtkm::cont::PartitionedDataSet CreateResult (const vtkm::cont::PartitionedDataSet &input, const vtkm::cont::PartitionedDataSet &resultPartitions) const
 Create the output data set for DoExecute. More...
 
template<typename FieldMapper >
vtkm::cont::PartitionedDataSet CreateResult (const vtkm::cont::PartitionedDataSet &input, const vtkm::cont::PartitionedDataSet &resultPartitions, FieldMapper &&fieldMapper) const
 Create the output data set for DoExecute. More...
 
template<typename FieldMapper >
vtkm::cont::DataSet CreateResult (const vtkm::cont::DataSet &inDataSet, const vtkm::cont::UnknownCellSet &resultCellSet, FieldMapper &&fieldMapper) const
 Create the output data set for DoExecute. More...
 
template<typename FieldMapper >
vtkm::cont::DataSet CreateResultCoordinateSystem (const vtkm::cont::DataSet &inDataSet, const vtkm::cont::UnknownCellSet &resultCellSet, const vtkm::cont::CoordinateSystem &resultCoordSystem, FieldMapper &&fieldMapper) const
 Create the output data set for DoExecute. More...
 
template<typename FieldMapper >
vtkm::cont::DataSet CreateResultCoordinateSystem (const vtkm::cont::DataSet &inDataSet, const vtkm::cont::UnknownCellSet &resultCellSet, const std::string &coordsName, const vtkm::cont::UnknownArrayHandle &coordsData, FieldMapper &&fieldMapper) const
 Create the output data set for DoExecute. More...
 
const vtkm::cont::FieldGetFieldFromDataSet (const vtkm::cont::DataSet &input) const
 Retrieve an input field from a vtkm::cont::DataSet object. More...
 
const vtkm::cont::FieldGetFieldFromDataSet (vtkm::IdComponent index, const vtkm::cont::DataSet &input) const
 Retrieve an input field from a vtkm::cont::DataSet object. More...
 
virtual vtkm::cont::DataSet DoExecute (const vtkm::cont::DataSet &inData)=0
 
virtual vtkm::cont::PartitionedDataSet DoExecutePartitions (const vtkm::cont::PartitionedDataSet &inData)
 
template<typename Functor , typename... Args>
void CastAndCallScalarField (const vtkm::cont::UnknownArrayHandle &fieldArray, Functor &&functor, Args &&... args) const
 Convenience method to get the array from a filter's input scalar field. More...
 
template<typename Functor , typename... Args>
void CastAndCallScalarField (const vtkm::cont::Field &field, Functor &&functor, Args &&... args) const
 Convenience method to get the array from a filter's input scalar field. More...
 
template<vtkm::IdComponent VecSize, typename Functor , typename... Args>
void CastAndCallVecField (const vtkm::cont::UnknownArrayHandle &fieldArray, Functor &&functor, Args &&... args) const
 Convenience method to get the array from a filter's input vector field. More...
 
template<vtkm::IdComponent VecSize, typename Functor , typename... Args>
void CastAndCallVecField (const vtkm::cont::Field &field, Functor &&functor, Args &&... args) const
 Convenience method to get the array from a filter's input vector field. More...
 
template<typename Functor , typename... Args>
void CastAndCallVariableVecField (const vtkm::cont::UnknownArrayHandle &fieldArray, Functor &&functor, Args &&... args) const
 This method is like CastAndCallVecField except that it can be used for a field of unknown vector size (or scalars). More...
 
template<typename Functor , typename... Args>
void CastAndCallVariableVecField (const vtkm::cont::Field &field, Functor &&functor, Args &&... args) const
 This method is like CastAndCallVecField except that it can be used for a field of unknown vector size (or scalars). More...
 

Protected Attributes

vtkm::cont::Invoker Invoke
 

Private Member Functions

template<typename FieldMapper >
void MapFieldsOntoOutput (const vtkm::cont::DataSet &input, const vtkm::filter::FieldSelection &fieldSelection, vtkm::cont::DataSet &output, FieldMapper &&fieldMapper) const
 
template<typename FieldMapper >
void MapFieldsOntoOutput (const vtkm::cont::PartitionedDataSet &input, const vtkm::filter::FieldSelection &fieldSelection, vtkm::cont::PartitionedDataSet &output, FieldMapper &&fieldMapper) const
 
virtual vtkm::Id DetermineNumberOfThreads (const vtkm::cont::PartitionedDataSet &input)
 
void ResizeIfNeeded (size_t index_st)
 

Private Attributes

vtkm::filter::FieldSelection FieldsToPass = vtkm::filter::FieldSelection::Mode::All
 
bool PassCoordinateSystems = true
 
bool RunFilterWithMultipleThreads = false
 
vtkm::Id NumThreadsPerCPU = 4
 
vtkm::Id NumThreadsPerGPU = 8
 
std::string OutputFieldName
 
std::vector< std::string > ActiveFieldNames
 
std::vector< vtkm::cont::Field::AssociationActiveFieldAssociation
 
std::vector< bool > UseCoordinateSystemAsField
 
std::vector< vtkm::IdActiveCoordinateSystemIndices
 

Detailed Description

Base class for all filters.

This is the base class for all filters. To add a new filter, one can subclass this and implement relevant methods.

FilterUsage Usage

To execute a filter, one typically calls the auto result = filter.Execute(input). Typical usage is as follows:

// create the concrete subclass (e.g. Contour).
// select fields to map to the output, if different from default which is to map all input
// fields.
contour.SetFieldsToPass({"var1", "var2"});
// execute the filter on vtkm::cont::DataSet.
vtkm::cont::DataSet dsInput = ...
auto outputDS = contour.Execute(dsInput);
// or, execute on a vtkm::cont::PartitionedDataSet
auto outputMB = contour.Execute(mbInput);

Execute methods take in the input DataSet or PartitionedDataSet to process and return the result. The type of the result is same as the input type, thus Execute(DataSet&) returns a DataSet while Execute(PartitionedDataSet&) returns a PartitionedDataSet.

Execute simply calls the pure virtual function DoExecute(DataSet&) which is the main extension point of the Filter interface. Filter developer needs to override DoExecute(DataSet) to implement the business logic of filtering operations on a single DataSet.

The default implementation of Execute(PartitionedDataSet&) is merely provided for convenience. Internally, it calls DoExecutePartitions(PartitionedDataSet) to iterate DataSets of a PartitionedDataSet and pass each individual DataSets to DoExecute(DataSet&), possibly in a multi-threaded setting. Developer of DoExecute(DataSet&) needs to indicate the thread-safeness of DoExecute(DataSet&) by overriding the CanThread() virtual method which by default returns true.

In the case that filtering on a PartitionedDataSet can not be simply implemented as a for-each loop on the component DataSets, filter implementor needs to override the DoExecutePartitions(PartitionedDataSet&). See the implementation of FilterParticleAdvection::Execute(PartitionedDataSet&) for an example.

Creating results and mapping fields

For subclasses that map input fields into output fields, the implementation of its DoExecute(DataSet&) should create the DataSet to be returned with a call to Filter::CreateResult or a similar method (such as Filter::CreateResultField).

VTKM_CONT DataSet SomeFilter::DoExecute(const vtkm::cont::DataSet& input)
{
outCellSet = ... // Generation of the new CellSet
// Mapper is a callable object (function object, lambda, etc.) that takes an input Field
// and maps it to an output Field and then add the output Field to the output DataSet
auto mapper = [](auto& outputDs, const auto& inputField) {
auto outputField = ... // Business logic for mapping input field to output field
output.AddField(outputField);
};
// This passes coordinate systems directly from input to output. If the points of
// the cell set change at all, they will have to be mapped by hand.
return this->CreateResult(input, outCellSet, mapper);
}

In addition to creating a new DataSet filled with the proper cell structure and coordinate systems, CreateResult iterates through each FieldToPass in the input DataSet and calls the FieldMapper to map the input Field to output Field. For simple filters that just pass on input fields to the output DataSet without any computation, an overload of CreateResult(const vtkm::cont::DataSet& input) is also provided as a convenience that uses the default mapper which trivially adds input Field to output DataSet (via a shallow copy).

FilterThreadSafety CanThread

By default, the implementation of DoExecute(DataSet&) should model a pure function, i.e. it does not have any mutable shared state. This makes it thread-safe by default and allows the default implementation of DoExecutePartitions(PartitionedDataSet&) to be simply a parallel for-each, thus facilitates multi-threaded execution without any lock.

Many legacy (VTKm 1.x) filter implementations needed to store states between the mesh generation phase and field mapping phase of filter execution, for example, parameters for field interpolation. The shared mutable states were mostly stored as mutable data members of the filter class (either in terms of ArrayHandle or some kind of Worket). The new filter interface, by combining the two phases into a single call to DoExecute(DataSet&), we have eliminated most of the cases that require such shared mutable states. New implementations of filters that require passing information between these two phases can now use local variables within the DoExecute(DataSet&). For example:

```cpp struct SharedState; // shared states between mesh generation and field mapping. VTKM_CONT DataSet ThreadSafeFilter::DoExecute(const vtkm::cont::DataSet& input) { // Mutable states that was a data member of the filter is now a local variable. // Each invocation of Execute(DataSet) in the multi-threaded execution of // Execute(PartitionedDataSet&) will have a copy of states on each thread's stack // thus making it thread-safe. SharedStates states;

vtkm::cont::CellSetExplicit<> cellSet; cellSet = ... // Generation of the new DataSet and store interpolation parameters in states

// Lambda capture of states, effectively passing the shared states to the Mapper. auto mapper = [&states](auto& outputDs, const auto& inputField) { auto outputField = ... // Use states for mapping input field to output field output.AddField(outputField); }; this->CreateOutput(input, cellSet, mapper);

return output; } ```

In the rare cases that filter implementation can not be made thread-safe, the implementation needs to override the CanThread() virtual method to return false. The default Execute(PartitionedDataSet&) implementation will fallback to a serial for loop execution.

FilterThreadScheduling DoExecute

The default multi-threaded execution of Execute(PartitionedDataSet&) uses a simple FIFO queue of DataSet and pool of worker threads. Implementation of Filter subclass can override the DoExecutePartitions(PartitionedDataSet) virtual method to provide implementation specific scheduling policy. The default number of worker threads in the pool are determined by the DetermineNumberOfThreads() virtual method using several backend dependent heuristic. Implementations of Filter subclass can also override DetermineNumberOfThreads() to provide implementation specific heuristic.

Constructor & Destructor Documentation

◆ Filter()

vtkm::filter::Filter::Filter ( )

◆ ~Filter()

virtual vtkm::filter::Filter::~Filter ( )
virtual

Member Function Documentation

◆ CanThread()

virtual bool vtkm::filter::Filter::CanThread ( ) const
virtual

Returns whether the filter can execute on partitions in concurrent threads.

If a derived class's implementation of DoExecute cannot run on multiple threads, then the derived class should override this method to return false.

Reimplemented in vtkm::filter::scalar_topology::ContourTreeUniformDistributed, vtkm::filter::scalar_topology::ContourTreeAugmented, vtkm::filter::entity_extraction::ExternalFaces, vtkm::filter::flow::FilterParticleAdvection, vtkm::filter::flow::LagrangianStructures, and vtkm::filter::flow::Lagrangian.

◆ CastAndCallScalarField() [1/2]

template<typename Functor , typename... Args>
void vtkm::filter::Filter::CastAndCallScalarField ( const vtkm::cont::Field field,
Functor &&  functor,
Args &&...  args 
) const
inlineprotected

Convenience method to get the array from a filter's input scalar field.

A field filter typically gets its input fields using the internal GetFieldFromDataSet. To use this field in a worklet, it eventually needs to be converted to an vtkm::cont::ArrayHandle. If the input field is limited to be a scalar field, then this method provides a convenient way to determine the correct array type. Like other CastAndCall methods, it takes as input a vtkm::cont::Field (or vtkm::cont::UnknownArrayHandle) and a function/functor to call with the appropriate vtkm::cont::ArrayHandle type.

◆ CastAndCallScalarField() [2/2]

template<typename Functor , typename... Args>
void vtkm::filter::Filter::CastAndCallScalarField ( const vtkm::cont::UnknownArrayHandle fieldArray,
Functor &&  functor,
Args &&...  args 
) const
inlineprotected

Convenience method to get the array from a filter's input scalar field.

A field filter typically gets its input fields using the internal GetFieldFromDataSet. To use this field in a worklet, it eventually needs to be converted to an vtkm::cont::ArrayHandle. If the input field is limited to be a scalar field, then this method provides a convenient way to determine the correct array type. Like other CastAndCall methods, it takes as input a vtkm::cont::Field (or vtkm::cont::UnknownArrayHandle) and a function/functor to call with the appropriate vtkm::cont::ArrayHandle type.

◆ CastAndCallVariableVecField() [1/2]

template<typename Functor , typename... Args>
void vtkm::filter::Filter::CastAndCallVariableVecField ( const vtkm::cont::Field field,
Functor &&  functor,
Args &&...  args 
) const
inlineprotected

This method is like CastAndCallVecField except that it can be used for a field of unknown vector size (or scalars).

This method will call the given functor with an vtkm::cont::ArrayHandleRecombineVec.

Note that there are limitations with using vtkm::cont::ArrayHandleRecombineVec within a worklet. Because the size of the vectors are not known at compile time, you cannot just create an intermediate vtkm::Vec of the correct size. Typically, you must allocate the output array (for example, with vtkm::cont::ArrayHandleRuntimeVec), and the worklet must iterate over the components and store them in the prealocated output.

◆ CastAndCallVariableVecField() [2/2]

template<typename Functor , typename... Args>
void vtkm::filter::Filter::CastAndCallVariableVecField ( const vtkm::cont::UnknownArrayHandle fieldArray,
Functor &&  functor,
Args &&...  args 
) const
inlineprotected

This method is like CastAndCallVecField except that it can be used for a field of unknown vector size (or scalars).

This method will call the given functor with an vtkm::cont::ArrayHandleRecombineVec.

Note that there are limitations with using vtkm::cont::ArrayHandleRecombineVec within a worklet. Because the size of the vectors are not known at compile time, you cannot just create an intermediate vtkm::Vec of the correct size. Typically, you must allocate the output array (for example, with vtkm::cont::ArrayHandleRuntimeVec), and the worklet must iterate over the components and store them in the prealocated output.

◆ CastAndCallVecField() [1/2]

template<vtkm::IdComponent VecSize, typename Functor , typename... Args>
void vtkm::filter::Filter::CastAndCallVecField ( const vtkm::cont::Field field,
Functor &&  functor,
Args &&...  args 
) const
inlineprotected

Convenience method to get the array from a filter's input vector field.

A field filter typically gets its input fields using the internal GetFieldFromDataSet. To use this field in a worklet, it eventually needs to be converted to an vtkm::cont::ArrayHandle. If the input field is limited to be a vector field with vectors of a specific size, then this method provides a convenient way to determine the correct array type. Like other CastAndCall methods, it takes as input a vtkm::cont::Field (or vtkm::cont::UnknownArrayHandle) and a function/functor to call with the appropriate vtkm::cont::ArrayHandle type. You also have to provide the vector size as the first template argument. For example CastAndCallVecField<3>(field, functor);.

◆ CastAndCallVecField() [2/2]

template<vtkm::IdComponent VecSize, typename Functor , typename... Args>
void vtkm::filter::Filter::CastAndCallVecField ( const vtkm::cont::UnknownArrayHandle fieldArray,
Functor &&  functor,
Args &&...  args 
) const
inlineprotected

Convenience method to get the array from a filter's input vector field.

A field filter typically gets its input fields using the internal GetFieldFromDataSet. To use this field in a worklet, it eventually needs to be converted to an vtkm::cont::ArrayHandle. If the input field is limited to be a vector field with vectors of a specific size, then this method provides a convenient way to determine the correct array type. Like other CastAndCall methods, it takes as input a vtkm::cont::Field (or vtkm::cont::UnknownArrayHandle) and a function/functor to call with the appropriate vtkm::cont::ArrayHandle type. You also have to provide the vector size as the first template argument. For example CastAndCallVecField<3>(field, functor);.

◆ CreateResult() [1/4]

vtkm::cont::DataSet vtkm::filter::Filter::CreateResult ( const vtkm::cont::DataSet inDataSet) const
protected

Create the output data set for DoExecute.

This form of CreateResult will create an output data set with the same cell structure and coordinate system as the input and pass all fields (as requested by the Filter state).

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with the cell set, coordinate system, and fields of inDataSet (as selected by the FieldsToPass state of the filter).

◆ CreateResult() [2/4]

template<typename FieldMapper >
vtkm::cont::DataSet vtkm::filter::Filter::CreateResult ( const vtkm::cont::DataSet inDataSet,
const vtkm::cont::UnknownCellSet resultCellSet,
FieldMapper &&  fieldMapper 
) const
inlineprotected

Create the output data set for DoExecute.

This form of CreateResult will create an output data set with the given CellSet. You must also provide a field mapper function, which is a function that takes the output DataSet being created and a Field from the input and then applies any necessary transformations to the field array and adds it to the DataSet.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultCellSetThe CellSet of the output will be set to this.
[in]fieldMapperA function or functor that takes a DataSet as its first argument and a Field as its second argument. The DataSet is the data being created and will eventually be returned by CreateResult. The Field comes from inDataSet. The function should map the Field to match resultCellSet and then add the resulting field to the DataSet. If the mapping is not possible, then the function should do nothing.

◆ CreateResult() [3/4]

vtkm::cont::PartitionedDataSet vtkm::filter::Filter::CreateResult ( const vtkm::cont::PartitionedDataSet input,
const vtkm::cont::PartitionedDataSet resultPartitions 
) const
protected

Create the output data set for DoExecute.

This form of CreateResult will create an output PartitionedDataSet with the same partitions and pass all PartitionedDataSet fields (as requested by the Filter state).

Parameters
[in]inputThe input data set being modified (usually the one passed into DoExecute).
[in]resultPartitionsThe output data created by the filter. Fields from the input are passed onto the return result partition as requested by the Filter state.

◆ CreateResult() [4/4]

template<typename FieldMapper >
vtkm::cont::PartitionedDataSet vtkm::filter::Filter::CreateResult ( const vtkm::cont::PartitionedDataSet input,
const vtkm::cont::PartitionedDataSet resultPartitions,
FieldMapper &&  fieldMapper 
) const
inlineprotected

Create the output data set for DoExecute.

This form of CreateResult will create an output PartitionedDataSet with the same partitions and pass all PartitionedDataSet fields (as requested by the Filter state).

Parameters
[in]inputThe input data set being modified (usually the one passed into DoExecute).
[in]resultPartitionsThe output data created by the filter. Fields from the input are passed onto the return result partition as requested by the Filter state.
[in]fieldMapperA function or functor that takes a PartitionedDataSet as its first argument and a Field as its second argument. The PartitionedDataSet is the data being created and will eventually be returned by CreateResult. The Field comes from input.

◆ CreateResultCoordinateSystem() [1/2]

template<typename FieldMapper >
vtkm::cont::DataSet vtkm::filter::Filter::CreateResultCoordinateSystem ( const vtkm::cont::DataSet inDataSet,
const vtkm::cont::UnknownCellSet resultCellSet,
const std::string &  coordsName,
const vtkm::cont::UnknownArrayHandle coordsData,
FieldMapper &&  fieldMapper 
) const
inlineprotected

Create the output data set for DoExecute.

This form of CreateResult will create an output data set with the given CellSet and CoordinateSystem. You must also provide a field mapper function, which is a function that takes the output DataSet being created and a Field from the input and then applies any necessary transformations to the field array and adds it to the DataSet.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultCellSetThe CellSet of the output will be set to this.
[in]coordsNameThe name of the coordinate system to be added to the output.
[in]coordsDataThe array containing the coordinates of the points.
[in]fieldMapperA function or functor that takes a DataSet as its first argument and a Field as its second argument. The DataSet is the data being created and will eventually be returned by CreateResult. The Field comes from inDataSet. The function should map the Field to match resultCellSet and then add the resulting field to the DataSet. If the mapping is not possible, then the function should do nothing.

◆ CreateResultCoordinateSystem() [2/2]

template<typename FieldMapper >
vtkm::cont::DataSet vtkm::filter::Filter::CreateResultCoordinateSystem ( const vtkm::cont::DataSet inDataSet,
const vtkm::cont::UnknownCellSet resultCellSet,
const vtkm::cont::CoordinateSystem resultCoordSystem,
FieldMapper &&  fieldMapper 
) const
inlineprotected

Create the output data set for DoExecute.

This form of CreateResult will create an output data set with the given CellSet and CoordinateSystem. You must also provide a field mapper function, which is a function that takes the output DataSet being created and a Field from the input and then applies any necessary transformations to the field array and adds it to the DataSet.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultCellSetThe CellSet of the output will be set to this.
[in]resultCoordSystemThis CoordinateSystem will be added to the output.
[in]fieldMapperA function or functor that takes a DataSet as its first argument and a Field as its second argument. The DataSet is the data being created and will eventually be returned by CreateResult. The Field comes from inDataSet. The function should map the Field to match resultCellSet and then add the resulting field to the DataSet. If the mapping is not possible, then the function should do nothing.

◆ CreateResultField() [1/2]

vtkm::cont::DataSet vtkm::filter::Filter::CreateResultField ( const vtkm::cont::DataSet inDataSet,
const std::string &  resultFieldName,
vtkm::cont::Field::Association  resultFieldAssociation,
const vtkm::cont::UnknownArrayHandle resultFieldArray 
) const
inlineprotected

Create the output data set for DoExecute

This form of CreateResult will create an output data set with the same cell structure and coordinate system as the input and pass all fields (as requested by the Filter state). Additionally, it will add a field matching the provided specifications to the result.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultFieldNameThe name of the field added to the returned DataSet.
[in]resultFieldAssociationThe association of the field (e.g. point or cell) added to the returned DataSet.
[in]resultFieldArrayAn array containing the data for the field added to the returned DataSet.

◆ CreateResultField() [2/2]

vtkm::cont::DataSet vtkm::filter::Filter::CreateResultField ( const vtkm::cont::DataSet inDataSet,
const vtkm::cont::Field resultField 
) const
protected

Create the output data set for DoExecute

This form of CreateResult will create an output data set with the same cell structure and coordinate system as the input and pass all fields (as requested by the Filter state). Additionally, it will add the provided field to the result.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultFieldA Field that is added to the returned DataSet.

◆ CreateResultFieldCell()

vtkm::cont::DataSet vtkm::filter::Filter::CreateResultFieldCell ( const vtkm::cont::DataSet inDataSet,
const std::string &  resultFieldName,
const vtkm::cont::UnknownArrayHandle resultFieldArray 
) const
inlineprotected

Create the output data set for DoExecute

This form of CreateResult will create an output data set with the same cell structure and coordinate system as the input and pass all fields (as requested by the Filter state). Additionally, it will add a cell field matching the provided specifications to the result.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultFieldNameThe name of the field added to the returned DataSet.
[in]resultFieldArrayAn array containing the data for the field added to the returned DataSet.

◆ CreateResultFieldPoint()

vtkm::cont::DataSet vtkm::filter::Filter::CreateResultFieldPoint ( const vtkm::cont::DataSet inDataSet,
const std::string &  resultFieldName,
const vtkm::cont::UnknownArrayHandle resultFieldArray 
) const
inlineprotected

Create the output data set for DoExecute

This form of CreateResult will create an output data set with the same cell structure and coordinate system as the input and pass all fields (as requested by the Filter state). Additionally, it will add a point field matching the provided specifications to the result.

Parameters
[in]inDataSetThe input data set being modified (usually the one passed into DoExecute). The returned DataSet is filled with fields of inDataSet (as selected by the FieldsToPass state of the filter).
[in]resultFieldNameThe name of the field added to the returned DataSet.
[in]resultFieldArrayAn array containing the data for the field added to the returned DataSet.

◆ DetermineNumberOfThreads()

virtual vtkm::Id vtkm::filter::Filter::DetermineNumberOfThreads ( const vtkm::cont::PartitionedDataSet input)
privatevirtual

◆ DoExecute()

virtual vtkm::cont::DataSet vtkm::filter::Filter::DoExecute ( const vtkm::cont::DataSet inData)
protectedpure virtual

Implemented in vtkm::filter::contour::AbstractContour, vtkm::filter::contour::Contour, vtkm::filter::contour::ContourMarchingCells, vtkm::filter::contour::ContourFlyingEdges, vtkm::filter::image_processing::ImageDifference, vtkm::filter::vector_analysis::SurfaceNormals, vtkm::filter::vector_analysis::Gradient, vtkm::filter::multi_block::MergeDataSets, vtkm::filter::scalar_topology::ContourTreeUniformDistributed, vtkm::filter::field_transform::PointTransform, vtkm::filter::field_transform::FieldToColors, vtkm::filter::scalar_topology::ContourTreeAugmented, vtkm::filter::field_transform::Warp, vtkm::filter::vector_analysis::DotProduct, vtkm::filter::entity_extraction::ExtractStructured, vtkm::filter::vector_analysis::CrossProduct, vtkm::filter::uncertainty::ContourUncertainUniform, vtkm::filter::mesh_info::MeshQuality, vtkm::filter::entity_extraction::Threshold, vtkm::filter::entity_extraction::ExtractGeometry, vtkm::filter::mesh_info::CellMeasures, vtkm::filter::scalar_topology::ContourTreeMesh3D, vtkm::filter::entity_extraction::GhostCellRemove, vtkm::filter::field_transform::GenerateIds, vtkm::filter::field_transform::LogValues, vtkm::filter::scalar_topology::ContourTreeMesh2D, vtkm::filter::entity_extraction::ExtractPoints, vtkm::filter::density_estimate::Histogram, vtkm::filter::resampling::HistSampling, vtkm::filter::field_transform::PointElevation, vtkm::filter::uncertainty::ContourUncertainUniformMonteCarlo, vtkm::filter::contour::MIRFilter, vtkm::filter::density_estimate::Statistics, vtkm::filter::entity_extraction::ExternalFaces, vtkm::filter::resampling::Probe, vtkm::filter::density_estimate::ParticleDensityNearestGridPoint, vtkm::filter::geometry_refinement::VertexClustering, vtkm::filter::density_estimate::ParticleDensityCloudInCell, vtkm::filter::geometry_refinement::SplitSharpEdges, vtkm::filter::mesh_info::MeshQualityScaledJacobian, vtkm::filter::geometry_refinement::ConvertToPointCloud, vtkm::filter::contour::ClipWithImplicitFunction, vtkm::filter::contour::ClipWithField, vtkm::filter::density_estimate::ContinuousScatterPlot, vtkm::filter::entity_extraction::ThresholdPoints, vtkm::filter::mesh_info::MeshQualityRelativeSizeSquared, vtkm::filter::mesh_info::MeshQualityWarpage, vtkm::filter::density_estimate::NDHistogram, vtkm::filter::field_transform::CompositeVectors, vtkm::filter::mesh_info::MeshQualityArea, vtkm::filter::mesh_info::MeshQualityVolume, vtkm::filter::mesh_info::MeshQualityShape, vtkm::filter::mesh_info::MeshQualityShapeAndSize, vtkm::filter::mesh_info::MeshQualitySkew, vtkm::filter::mesh_info::MeshQualityTaper, vtkm::filter::geometry_refinement::Tube, vtkm::filter::mesh_info::MeshQualityMaxAngle, vtkm::filter::mesh_info::MeshQualityMinAngle, vtkm::filter::mesh_info::MeshQualityStretch, vtkm::filter::contour::Slice, vtkm::filter::geometry_refinement::Shrink, vtkm::filter::mesh_info::MeshQualityAspectGamma, vtkm::filter::mesh_info::MeshQualityAspectRatio, vtkm::filter::mesh_info::MeshQualityShear, vtkm::filter::mesh_info::MeshQualityCondition, vtkm::filter::mesh_info::MeshQualityDiagonalRatio, vtkm::filter::mesh_info::MeshQualityOddy, vtkm::filter::contour::SliceMultiple, vtkm::filter::entity_extraction::Mask, vtkm::filter::field_transform::CylindricalCoordinateTransform, vtkm::filter::density_estimate::Entropy, vtkm::filter::entity_extraction::MaskPoints, vtkm::filter::field_transform::SphericalCoordinateTransform, vtkm::filter::connected_components::ImageConnectivity, vtkm::filter::image_processing::ImageMedian, vtkm::filter::mesh_info::MeshQualityDimension, vtkm::filter::mesh_info::MeshQualityJacobian, vtkm::filter::mesh_info::MeshQualityMaxDiagonal, vtkm::filter::mesh_info::MeshQualityMinDiagonal, vtkm::filter::connected_components::CellSetConnectivity, vtkm::filter::zfp::ZFPCompressor2D, vtkm::filter::zfp::ZFPCompressor3D, vtkm::filter::zfp::ZFPDecompressor1D, vtkm::filter::zfp::ZFPDecompressor2D, vtkm::filter::zfp::ZFPDecompressor3D, vtkm::filter::field_conversion::CellAverage, vtkm::filter::field_conversion::PointAverage, vtkm::filter::vector_analysis::VectorMagnitude, vtkm::filter::geometry_refinement::Tetrahedralize, vtkm::filter::geometry_refinement::Triangulate, vtkm::filter::image_processing::ComputeMoments, vtkm::filter::density_estimate::NDEntropy, vtkm::filter::zfp::ZFPCompressor1D, vtkm::filter::flow::FilterParticleAdvection, vtkm::filter::flow::LagrangianStructures, vtkm::filter::flow::Lagrangian, vtkm::filter::clean_grid::CleanGrid, vtkm::filter::flow::StreamSurface, vtkm::filter::mesh_info::GhostCellClassify, vtkm::filter::scalar_topology::DistributedBranchDecompositionFilter, vtkm::filter::scalar_topology::SelectTopVolumeContoursFilter, and vtkm::filter::multi_block::AmrArrays.

◆ DoExecutePartitions()

virtual vtkm::cont::PartitionedDataSet vtkm::filter::Filter::DoExecutePartitions ( const vtkm::cont::PartitionedDataSet inData)
protectedvirtual

◆ Execute() [1/2]

vtkm::cont::DataSet vtkm::filter::Filter::Execute ( const vtkm::cont::DataSet input)

Executes the filter on the input and produces a result dataset.

On success, this the dataset produced. On error, vtkm::cont::ErrorExecution will be thrown.

◆ Execute() [2/2]

vtkm::cont::PartitionedDataSet vtkm::filter::Filter::Execute ( const vtkm::cont::PartitionedDataSet input)

Executes the filter on the input PartitionedDataSet and produces a result PartitionedDataSet.

On success, this the dataset produced. On error, vtkm::cont::ErrorExecution will be thrown.

◆ GetActiveCoordinateSystemIndex()

vtkm::Id vtkm::filter::Filter::GetActiveCoordinateSystemIndex ( vtkm::IdComponent  index = 0) const
inline

Specifies the coordinate system index to make active to use when processing the input vtkm::cont::DataSet.

This is used primarily by the Filter to select the coordinate system to use as a field when UseCoordinateSystemAsField is true.

◆ GetActiveFieldAssociation()

vtkm::cont::Field::Association vtkm::filter::Filter::GetActiveFieldAssociation ( vtkm::IdComponent  index = 0) const
inline

Specifies a field to operate on.

The number of input fields (or whether the filter operates on input fields at all) is specific to each particular filter.

◆ GetActiveFieldName()

const std::string& vtkm::filter::Filter::GetActiveFieldName ( vtkm::IdComponent  index = 0) const
inline

Specifies a field to operate on.

The number of input fields (or whether the filter operates on input fields at all) is specific to each particular filter.

◆ GetFieldFromDataSet() [1/2]

const vtkm::cont::Field& vtkm::filter::Filter::GetFieldFromDataSet ( const vtkm::cont::DataSet input) const
inlineprotected

Retrieve an input field from a vtkm::cont::DataSet object.

When a filter operates on fields, it should use this method to get the input fields that the use has selected with SetActiveField() and related methods.

◆ GetFieldFromDataSet() [2/2]

const vtkm::cont::Field& vtkm::filter::Filter::GetFieldFromDataSet ( vtkm::IdComponent  index,
const vtkm::cont::DataSet input 
) const
inlineprotected

Retrieve an input field from a vtkm::cont::DataSet object.

When a filter operates on fields, it should use this method to get the input fields that the use has selected with SetActiveField() and related methods.

◆ GetFieldsToPass() [1/2]

vtkm::filter::FieldSelection& vtkm::filter::Filter::GetFieldsToPass ( )
inline

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ GetFieldsToPass() [2/2]

const vtkm::filter::FieldSelection& vtkm::filter::Filter::GetFieldsToPass ( ) const
inline

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ GetNumberOfActiveFields()

vtkm::IdComponent vtkm::filter::Filter::GetNumberOfActiveFields ( ) const
inline

Return the number of active fields currently set.

The general interface to Filter allows a user to set an arbitrary number of active fields (indexed 0 and on). This method returns the number of active fields that are set. Note that the filter implementation is free to ignore any active fields it does not support. Also note that an active field can be set to be either a named field or a coordinate system.

◆ GetOutputFieldName()

const std::string& vtkm::filter::Filter::GetOutputFieldName ( ) const
inline

Specifies the name of the output field generated.

Not all filters create an output field.

◆ GetPassCoordinateSystems()

bool vtkm::filter::Filter::GetPassCoordinateSystems ( ) const
inline

Specify whether to always pass coordinate systems.

vtkm::cont::CoordinateSystems in a DataSet are really just point fields marked as being a coordinate system. Thus, a coordinate system is passed if and only if the associated field is passed.

By default, the filter will pass all fields associated with a coordinate system regardless of the FieldsToPass marks the field as passing. If this option is set to false, then coordinate systems will only be passed if it is marked so by FieldsToPass.

◆ GetRunMultiThreadedFilter()

bool vtkm::filter::Filter::GetRunMultiThreadedFilter ( ) const
inline

◆ GetThreadsPerCPU()

vtkm::Id vtkm::filter::Filter::GetThreadsPerCPU ( ) const
inline

◆ GetThreadsPerGPU()

vtkm::Id vtkm::filter::Filter::GetThreadsPerGPU ( ) const
inline

◆ GetUseCoordinateSystemAsField()

bool vtkm::filter::Filter::GetUseCoordinateSystemAsField ( vtkm::IdComponent  index = 0) const
inline

Specifies whether to use point coordinates as the input field.

When true, the values for the active field are ignored and the active coordinate system is used instead.

◆ MapFieldsOntoOutput() [1/2]

template<typename FieldMapper >
void vtkm::filter::Filter::MapFieldsOntoOutput ( const vtkm::cont::DataSet input,
const vtkm::filter::FieldSelection fieldSelection,
vtkm::cont::DataSet output,
FieldMapper &&  fieldMapper 
) const
inlineprivate

◆ MapFieldsOntoOutput() [2/2]

template<typename FieldMapper >
void vtkm::filter::Filter::MapFieldsOntoOutput ( const vtkm::cont::PartitionedDataSet input,
const vtkm::filter::FieldSelection fieldSelection,
vtkm::cont::PartitionedDataSet output,
FieldMapper &&  fieldMapper 
) const
inlineprivate

◆ ResizeIfNeeded()

void vtkm::filter::Filter::ResizeIfNeeded ( size_t  index_st)
private

◆ SetActiveCoordinateSystem() [1/2]

void vtkm::filter::Filter::SetActiveCoordinateSystem ( vtkm::Id  coord_idx)
inline

Specifies the coordinate system index to make active to use when processing the input vtkm::cont::DataSet.

This is used primarily by the Filter to select the coordinate system to use as a field when UseCoordinateSystemAsField is true.

◆ SetActiveCoordinateSystem() [2/2]

void vtkm::filter::Filter::SetActiveCoordinateSystem ( vtkm::IdComponent  index,
vtkm::Id  coord_idx 
)
inline

Specifies the coordinate system index to make active to use when processing the input vtkm::cont::DataSet.

This is used primarily by the Filter to select the coordinate system to use as a field when UseCoordinateSystemAsField is true.

◆ SetActiveField() [1/2]

void vtkm::filter::Filter::SetActiveField ( const std::string &  name,
vtkm::cont::Field::Association  association = vtkm::cont::Field::Association::Any 
)
inline

Specifies a field to operate on.

The number of input fields (or whether the filter operates on input fields at all) is specific to each particular filter.

◆ SetActiveField() [2/2]

void vtkm::filter::Filter::SetActiveField ( vtkm::IdComponent  index,
const std::string &  name,
vtkm::cont::Field::Association  association = vtkm::cont::Field::Association::Any 
)
inline

Specifies a field to operate on.

The number of input fields (or whether the filter operates on input fields at all) is specific to each particular filter.

◆ SetFieldsToPass() [1/7]

void vtkm::filter::Filter::SetFieldsToPass ( const std::string &  fieldname,
vtkm::cont::Field::Association  association,
vtkm::filter::FieldSelection::Mode  mode = vtkm::filter::FieldSelection::Mode::Select 
)

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ SetFieldsToPass() [2/7]

void vtkm::filter::Filter::SetFieldsToPass ( const std::string &  fieldname,
vtkm::filter::FieldSelection::Mode  mode 
)
inline

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ SetFieldsToPass() [3/7]

void vtkm::filter::Filter::SetFieldsToPass ( const vtkm::filter::FieldSelection fieldsToPass)

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ SetFieldsToPass() [4/7]

void vtkm::filter::Filter::SetFieldsToPass ( const vtkm::filter::FieldSelection fieldsToPass,
vtkm::filter::FieldSelection::Mode  mode 
)

◆ SetFieldsToPass() [5/7]

void vtkm::filter::Filter::SetFieldsToPass ( std::initializer_list< std::pair< std::string, vtkm::cont::Field::Association >>  fields,
vtkm::filter::FieldSelection::Mode  mode = vtkm::filter::FieldSelection::Mode::Select 
)

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ SetFieldsToPass() [6/7]

void vtkm::filter::Filter::SetFieldsToPass ( std::initializer_list< std::string >  fields,
vtkm::filter::FieldSelection::Mode  mode = vtkm::filter::FieldSelection::Mode::Select 
)

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ SetFieldsToPass() [7/7]

void vtkm::filter::Filter::SetFieldsToPass ( vtkm::filter::FieldSelection &&  fieldsToPass)

Specify which fields get passed from input to output.

After a filter successfully executes and returns a new data set, fields are mapped from input to output. Depending on what operation the filter does, this could be a simple shallow copy of an array, or it could be a computed operation. You can control which fields are passed (and equivalently which are not) with this parameter.

By default, all fields are passed during execution.

◆ SetInvoker()

void vtkm::filter::Filter::SetInvoker ( vtkm::cont::Invoker  inv)
inline

Specify the vtkm::cont::Invoker to be used to execute worklets by this filter instance.

Overriding the default allows callers to control which device adapters a filter uses.

◆ SetOutputFieldName()

void vtkm::filter::Filter::SetOutputFieldName ( const std::string &  name)
inline

Specifies the name of the output field generated.

Not all filters create an output field.

◆ SetPassCoordinateSystems()

void vtkm::filter::Filter::SetPassCoordinateSystems ( bool  flag)
inline

Specify whether to always pass coordinate systems.

vtkm::cont::CoordinateSystems in a DataSet are really just point fields marked as being a coordinate system. Thus, a coordinate system is passed if and only if the associated field is passed.

By default, the filter will pass all fields associated with a coordinate system regardless of the FieldsToPass marks the field as passing. If this option is set to false, then coordinate systems will only be passed if it is marked so by FieldsToPass.

◆ SetRunMultiThreadedFilter()

void vtkm::filter::Filter::SetRunMultiThreadedFilter ( bool  val)
inline

◆ SetThreadsPerCPU()

void vtkm::filter::Filter::SetThreadsPerCPU ( vtkm::Id  numThreads)
inline

◆ SetThreadsPerGPU()

void vtkm::filter::Filter::SetThreadsPerGPU ( vtkm::Id  numThreads)
inline

◆ SetUseCoordinateSystemAsField() [1/2]

void vtkm::filter::Filter::SetUseCoordinateSystemAsField ( bool  val)
inline

Specifies whether to use point coordinates as the input field.

When true, the values for the active field are ignored and the active coordinate system is used instead.

◆ SetUseCoordinateSystemAsField() [2/2]

void vtkm::filter::Filter::SetUseCoordinateSystemAsField ( vtkm::IdComponent  index,
bool  val 
)
inline

Specifies whether to use point coordinates as the input field.

When true, the values for the active field are ignored and the active coordinate system is used instead.

Member Data Documentation

◆ ActiveCoordinateSystemIndices

std::vector<vtkm::Id> vtkm::filter::Filter::ActiveCoordinateSystemIndices
private

◆ ActiveFieldAssociation

std::vector<vtkm::cont::Field::Association> vtkm::filter::Filter::ActiveFieldAssociation
private

◆ ActiveFieldNames

std::vector<std::string> vtkm::filter::Filter::ActiveFieldNames
private

◆ FieldsToPass

vtkm::filter::FieldSelection vtkm::filter::Filter::FieldsToPass = vtkm::filter::FieldSelection::Mode::All
private

◆ Invoke

vtkm::cont::Invoker vtkm::filter::Filter::Invoke
protected

◆ NumThreadsPerCPU

vtkm::Id vtkm::filter::Filter::NumThreadsPerCPU = 4
private

◆ NumThreadsPerGPU

vtkm::Id vtkm::filter::Filter::NumThreadsPerGPU = 8
private

◆ OutputFieldName

std::string vtkm::filter::Filter::OutputFieldName
private

◆ PassCoordinateSystems

bool vtkm::filter::Filter::PassCoordinateSystems = true
private

◆ RunFilterWithMultipleThreads

bool vtkm::filter::Filter::RunFilterWithMultipleThreads = false
private

◆ UseCoordinateSystemAsField

std::vector<bool> vtkm::filter::Filter::UseCoordinateSystemAsField
private

The documentation for this class was generated from the following file:
vtkm::cont::DataSet
Contains and manages the geometric data structures that VTK-m operates on.
Definition: DataSet.h:57
vtkm::cont::UnknownCellSet
A CellSet of an unknown type.
Definition: UnknownCellSet.h:48
vtkm::filter::Filter::SetFieldsToPass
void SetFieldsToPass(const vtkm::filter::FieldSelection &fieldsToPass)
Specify which fields get passed from input to output.
VTKM_CONT
#define VTKM_CONT
Definition: ExportMacros.h:57
vtkm::filter::contour::Contour
Generate contours or isosurfaces from a region of space.
Definition: Contour.h:35
vtkm::filter::Filter::Execute
vtkm::cont::DataSet Execute(const vtkm::cont::DataSet &input)
Executes the filter on the input and produces a result dataset.
vtkm::cont::PartitionedDataSet
Comprises a set of vtkm::cont::DataSet objects.
Definition: PartitionedDataSet.h:26
vtkm::filter::Filter::CreateResult
vtkm::cont::DataSet CreateResult(const vtkm::cont::DataSet &inDataSet) const
Create the output data set for DoExecute.