2.7. Provided Filters

VTK‑m comes with the implementation of many filters. Filters in VTK‑m are divided into a collection of modules, each with its own namespace and library. This section is organized by each filter module, each of which contains one or more filters that are related to each other.

Note that this is not an exhaustive list of filters available in VTK‑m. More can be found in the namespaces under vtkm::filter (and likewise the subdirectories under vtkm/filter in the VTK‑m source.

2.7.1. Cleaning Grids

The vtkm::filter::clean_grid module contains filters that resolve issues with mesh structure. This could include finding and merging coincident points, removing degenerate cells, or converting the grid to a known type.

2.7.1.1. Clean Grid

vtkm::filter::clean_grid::CleanGrid is a filter that converts a cell set to an explicit representation and potentially removes redundant or unused data. It does this by iterating over all cells in the data set, and for each one creating the explicit cell representation that is stored in the output. (Explicit cell sets are described in Section 2.4.2.2 (Explicit Cell Sets).) One benefit of using vtkm::filter::clean_grid::CleanGrid is that it can optionally remove unused points and combine coincident points. Another benefit is that the resulting cell set will be of a known specific type.

Common Errors

The result of vtkm::filter::clean_grid::CleanGrid is not necessarily smaller, memory-wise, than its input. For example, “cleaning” a data set with a structured topology will actually result in a data set that requires much more memory to store an explicit topology.

class CleanGrid : public vtkm::filter::Filter

Clean a mesh to an unstructured grid.

This filter converts the cells of its input to an explicit representation and potentially removes redundant or unused data. The newly constructed data set will have the same cells as the input and the topology will be stored in a vtkm::cont::CellSetExplicit<>. The filter will also optionally remove all unused points.

Note that the result of CleanGrid is not necessarily smaller than the input. For example, “cleaning” a data set with a vtkm::cont::CellSetStructured topology will actually result in a much larger data set.

CleanGrid can optionally merge close points. The closeness of points is determined by the coordinate system. If there are multiple coordinate systems, the desired coordinate system can be selected with the SetActiveCoordinateSystem() method.

Public Functions

inline bool GetCompactPointFields() const

When the CompactPointFields flag is true, the filter will identify and remove any points that are not used by the topology.

This is on by default.

inline void SetCompactPointFields(bool flag)

When the CompactPointFields flag is true, the filter will identify and remove any points that are not used by the topology.

This is on by default.

inline bool GetMergePoints() const

When the MergePoints flag is true, the filter will identify any coincident points and merge them together.

The distance two points can be to considered coincident is set with the tolerance flags. This is on by default.

inline void SetMergePoints(bool flag)

When the MergePoints flag is true, the filter will identify any coincident points and merge them together.

The distance two points can be to considered coincident is set with the tolerance flags. This is on by default.

inline vtkm::Float64 GetTolerance() const

Defines the tolerance used when determining whether two points are considered coincident.

Because floating point parameters have limited precision, point coordinates that are essentially the same might not be bit-wise exactly the same. Thus, the CleanGrid filter has the ability to find and merge points that are close but perhaps not exact. If the ToleranceIsAbsolute flag is false (the default), then this tolerance is scaled by the diagonal of the points.

inline void SetTolerance(vtkm::Float64 tolerance)

Defines the tolerance used when determining whether two points are considered coincident.

Because floating point parameters have limited precision, point coordinates that are essentially the same might not be bit-wise exactly the same. Thus, the CleanGrid filter has the ability to find and merge points that are close but perhaps not exact. If the ToleranceIsAbsolute flag is false (the default), then this tolerance is scaled by the diagonal of the points.

inline bool GetToleranceIsAbsolute() const

When ToleranceIsAbsolute is false (the default) then the tolerance is scaled by the diagonal of the bounds of the dataset.

If true, then the tolerance is taken as the actual distance to use.

inline void SetToleranceIsAbsolute(bool flag)

When ToleranceIsAbsolute is false (the default) then the tolerance is scaled by the diagonal of the bounds of the dataset.

If true, then the tolerance is taken as the actual distance to use.

inline bool GetRemoveDegenerateCells() const

When RemoveDegenerateCells is true (the default), then CleanGrid will look for repeated points in cells and, if the repeated points cause the cell to drop dimensionality, the cell is removed.

This is particularly useful when point merging is on as this operation can create degenerate cells.

inline void SetRemoveDegenerateCells(bool flag)

When RemoveDegenerateCells is true (the default), then CleanGrid will look for repeated points in cells and, if the repeated points cause the cell to drop dimensionality, the cell is removed.

This is particularly useful when point merging is on as this operation can create degenerate cells.

inline bool GetFastMerge() const

When FastMerge is true (the default), some corners are cut when computing coincident points.

The point merge will go faster but the tolerance will not be strictly followed.

inline void SetFastMerge(bool flag)

When FastMerge is true (the default), some corners are cut when computing coincident points.

The point merge will go faster but the tolerance will not be strictly followed.

2.7.2. Connected Components

Connected components in a mesh are groups of mesh elements that are connected together in some way. For example, if two cells are neighbors, then they are in the same component. Likewise, a cell is also in the same component as its neighbor’s neighbors as well as their neighbors and so on. Connected components help identify when features in a simulation fragment or meld.

The vtkm::filter::connected_components module contains filters that find groups of cells that are connected. There are different ways to define what it means to be connected. One way is to use the topological connections of the cells. That is, two cells that share a point, edge, or face are connected. Another way is to use a field that classifies each cell, and cells are only connected if they have the same classification.

2.7.2.1. Cell Connectivity

The vtkm::filter::connected_components::CellSetConnectivity filter finds groups of cells that are connected together through their topology.

class CellSetConnectivity : public vtkm::filter::Filter

Finds and labels groups of cells that are connected together through their topology.

Two cells are considered connected if they share an edge. CellSetConnectivity identifies some number of components and assigns each component a unique integer.

The result of the filter is a cell field of type vtkm::Id with the default name of “component” (which can be changed with the SetOutputFieldName method). Each entry in the cell field will be a number that identifies to which component the cell belongs.

2.7.2.2. Classification Field on Image Data

The vtkm::filter::connected_components::ImageConnectivity filter finds groups of points that have the same field value and are connected together through their topology.

class ImageConnectivity : public vtkm::filter::Filter

2.7.3. Contouring

The vtkm::filter::contour module contains filters that extract regions that match some field or spatial criteria. Unlike entity extraction filters (Section 2.7.5), the geometry will be clipped or sliced to extract the exact matching region. (In contrast, entity extraction filters will pull unmodified points, edges, faces, or cells from the input.)

2.7.3.1. Contour

Contouring is one of the most fundamental filters in scientific visualization. A contour is the locus where a field is equal to a particular value. A topographic map showing curves of various elevations often used when hiking in hilly regions is an example of contours of an elevation field in 2 dimensions. Extended to 3 dimensions, a contour gives a surface. Thus, a contour is often called an isosurface. The contouring/isosurface algorithm is implemented by vtkm::filter::contour::Contour.

class Contour : public vtkm::filter::contour::AbstractContour

Generate contours or isosurfaces from a region of space.

Contour takes as input a mesh, often a volume, and generates on output one or more surfaces where a field equals a specified value.

This filter implements multiple algorithms for contouring, and the best algorithm will be selected based on the type of the input.

The scalar field to extract the contour from is selected with the SetActiveField() and related methods.

Subclassed by vtkm::filter::contour::Slice, vtkm::filter::contour::SliceMultiple

vtkm::filter::contour::Contour also inherits the following methods.

inline void vtkm::filter::contour::AbstractContour::SetIsoValue(vtkm::Float64 v)

Set a field value on which to extract a contour.

This form of the method is usually used when only one contour is being extracted.

inline void vtkm::filter::contour::AbstractContour::SetIsoValue(vtkm::Id index, vtkm::Float64 v)

Set a field value on which to extract a contour.

This form is used to specify multiple contours. The method is called multiple times with different index parameters.

inline void vtkm::filter::contour::AbstractContour::SetIsoValues(const std::vector<vtkm::Float64> &values)

Set multiple iso values at once.

The iso values can be specified as either a std::vector or an initializer list. So, both

std::vector<vtkm::Float64> isovalues = { 0.2, 0.5, 0.7 };
contour.SetIsoValues(isovalues);

and

contour.SetIsoValues({ 0.2, 0.5, 0.7 });

work.

inline vtkm::Float64 vtkm::filter::contour::AbstractContour::GetIsoValue(vtkm::Id index = 0) const

Return a value used to contour the mesh.

inline void vtkm::filter::contour::AbstractContour::SetGenerateNormals(bool flag)

Set whether normals should be generated.

Normals are used in shading calculations during rendering and can make the surface appear more smooth.

Off by default.

inline bool vtkm::filter::contour::AbstractContour::GetGenerateNormals() const

Get whether normals should be generated.

inline void vtkm::filter::contour::AbstractContour::SetComputeFastNormals(bool flag)

Set whether the fast path should be used for normals computation.

When this flag is off (the default), the generated normals are based on the gradient of the field being contoured and can be quite expensive to compute. When the flag is on, a faster method that computes the normals based on the faces of the isosurface mesh is used, but the normals do not look as good as the gradient based normals.

This flag has no effect if SetGenerateNormals is false.

inline bool vtkm::filter::contour::AbstractContour::GetComputeFastNormals() const

Get whether the fast path should be used for normals computation.

inline void vtkm::filter::contour::AbstractContour::SetNormalArrayName(const std::string &name)

Set the name of the field for the generated normals.

inline const std::string &vtkm::filter::contour::AbstractContour::GetNormalArrayName() const

Get the name of the field for the generated normals.

inline void vtkm::filter::contour::AbstractContour::SetMergeDuplicatePoints(bool on)

Set whether the points generated should be unique for every triangle or will duplicate points be merged together.

Duplicate points are identified by the unique edge it was generated from.

Because the contour filter (like all filters in VTK-m) runs in parallel, parallel threads can (and often do) create duplicate versions of points. When this flag is set to true, a secondary operation will find all duplicated points and combine them together. If false, points will be duplicated. In addition to requiring more storage, duplicated points mean that triangles next to each other will not be considered adjecent to subsequent filters.

inline bool vtkm::filter::contour::AbstractContour::GetMergeDuplicatePoints()

Get whether the points generated should be unique for every triangle or will duplicate points be merged together.

1  vtkm::filter::contour::Contour contour;
2
3  contour.SetActiveField("pointvar");
4  contour.SetIsoValue(10.0);
5
6  vtkm::cont::DataSet isosurface = contour.Execute(inData);

2.7.3.2. Slice

A slice operation intersects a mesh with a surface. The vtkm::filter::contour::Slice filter uses a vtkm::ImplicitFunctionGeneral to specify an implicit surface to slice on. A plane is a common thing to slice on, but other surfaces are available. See Chapter 2.12 (Implicit Functions) for information on implicit functions.

class Slice : public vtkm::filter::contour::Contour

Intersect a mesh with an implicit surface.

This filter accepts a vtkm::ImplicitFunction that defines the surface to slice on. A vtkm::Plane is a common function to use that cuts the mesh along a plane.

Public Functions

inline void SetImplicitFunction(const vtkm::ImplicitFunctionGeneral &func)

Set the implicit function that is used to perform the slicing.

Only a limited number of implicit functions are supported. See vtkm::ImplicitFunctionGeneral for information on which ones.

inline const vtkm::ImplicitFunctionGeneral &GetImplicitFunction() const

Get the implicit function that us used to perform the slicing.

The vtkm::filter::contour::Slice filter inherits from the vtkm::filter::contour::Contour, uses its implementation to extract the slices, and several of the inherited methods are useful including vtkm::filter::contour::AbstractContour::SetGenerateNormals(), vtkm::filter::contour::AbstractContour::GetGenerateNormals(), vtkm::filter::contour::AbstractContour::SetComputeFastNormals(), vtkm::filter::contour::AbstractContour::GetComputeFastNormals(), vtkm::filter::contour::AbstractContour::SetNormalArrayName(), vtkm::filter::contour::AbstractContour::GetNormalArrayName(), vtkm::filter::contour::AbstractContour::SetMergeDuplicatePoints(), vtkm::filter::contour::AbstractContour::GetMergeDuplicatePoints(), vtkm::filter::Field::SetActiveCoordinateSystem(), and vtkm::filter::Field::GetActiveCoordinateSystemIndex().

2.7.3.3. Clip with Field

Clipping is an operation that removes regions from the data set based on a user-provided value or function. The vtkm::filter::contour::ClipWithField filter takes a clip value as an argument and removes regions where a named scalar field is below (or above) that value. (A companion filter that discards a region of the data based on an implicit function is described later.)

The result of vtkm::filter::contour::ClipWithField is a volume. If a cell has field values at its vertices that are all below the specified value, then it will be discarded entirely. Likewise, if a cell has field values at its vertices that are all above the specified value, then it will be retained in its entirety. If a cell has some vertices with field values below the specified value and some above, then the cell will be split into the portions above the value (which will be retained) and the portions below the value (which will be discarded).

This operation is sometimes called an isovolume because it extracts the volume of a mesh that is inside the iso-region of a scalar. This is in contrast to an isosurface, which extracts only the surface of that iso-value. That said, a more appropriate name is interval volume as the volume is defined by a range of values, not a single “iso” value.

vtkm::filter::contour::ClipWithField is also similar to a threshold operation, which extracts cells based on the value of field. The difference is that threshold will either keep or remove entire cells based on the field values whereas clip with carve cells that straddle the valid regions. See Section 2.7.5.6 (Threshold) for information on threshold extraction.

class ClipWithField : public vtkm::filter::Filter

Clip a dataset using a field.

Clip a dataset using a given field value. All points that are less than that value are considered outside, and will be discarded. All points that are greater are kept.

To select the scalar field, use the SetActiveField() and related methods.

Public Functions

inline void SetClipValue(vtkm::Float64 value)

Specifies the field value for the clip operation.

Regions where the active field is less than this value are clipped away from each input cell.

inline void SetInvertClip(bool invert)

Specifies if the result for the clip filter should be inverted.

If set to false (the default), regions where the active field is less than the specified clip value are removed. If set to true, regions where the active field is more than the specified clip value are removed.

inline vtkm::Float64 GetClipValue() const

Specifies the field value for the clip operation.

inline bool GetInvertClip() const

Specifies if the result for the clip filter should be inverted.

1  // Create an instance of a clip filter that discards all regions with scalar
2  // value less than 25.
3  vtkm::filter::contour::ClipWithField clip;
4  clip.SetClipValue(25.0);
5  clip.SetActiveField("pointvar");
6
7  // Execute the clip filter
8  vtkm::cont::DataSet outData = clip.Execute(inData);

2.7.3.4. Clip with Implicit Function

The vtkm::filter::contour::ClipWithImplicitFunction function takes an implicit function and removes all parts of the data that are inside (or outside) that function. See Chapter 2.12 (Implicit Functions) for more detail on how implicit functions are represented in VTK‑m. A companion filter that discards a region of the data based on the value of a scalar field is described in Section 2.7.5.2 (Extract Geometry).

The result of vtkm::filter::contour::ClipWithImplicitFunction is a volume. If a cell has its vertices positioned all outside the implicit function, then it will be discarded entirely. Likewise, if a cell its vertices all inside the implicit function, then it will be retained in its entirety. If a cell has some vertices inside the implicit function and some outside, then the cell will be split into the portions inside (which will be retained) and the portions outside (which will be discarded).

class ClipWithImplicitFunction : public vtkm::filter::Filter

Clip a dataset using an implicit function.

Clip a dataset using a given implicit function value, such as vtkm::Sphere or vtkm::Frustum. The implicit function uses the point coordinates as its values. If there is more than one coordinate system in the input vtkm::cont::DataSet, it can be selected with SetActiveCoordinateSystem().

Public Functions

inline void SetImplicitFunction(const vtkm::ImplicitFunctionGeneral &func)

Specifies the implicit function to be used to perform the clip operation.

Only a limited number of implicit functions are supported. See vtkm::ImplicitFunctionGeneral for information on which ones.

inline void SetInvertClip(bool invert)

Specifies whether the result of the clip filter should be inverted.

If set to false (the default), all regions where the implicit function is negative will be removed. If set to true, all regions where the implicit function is positive will be removed.

inline const vtkm::ImplicitFunctionGeneral &GetImplicitFunction() const

Specifies the implicit function to be used to perform the clip operation.

In the example provided below the vtkm::Sphere implicit function is used. This function evaluates to a negative value if points from the original dataset occur within the sphere, evaluates to 0 if the points occur on the surface of the sphere, and evaluates to a positive value if the points occur outside the sphere.

 1  // Parameters needed for implicit function
 2  vtkm::Sphere implicitFunction(vtkm::make_Vec(1, 0, 1), 0.5);
 3
 4  // Create an instance of a clip filter with this implicit function.
 5  vtkm::filter::contour::ClipWithImplicitFunction clip;
 6  clip.SetImplicitFunction(implicitFunction);
 7
 8  // By default, ClipWithImplicitFunction will remove everything inside the sphere.
 9  // Set the invert clip flag to keep the inside of the sphere and remove everything
10  // else.
11  clip.SetInvertClip(true);
12
13  // Execute the clip filter
14  vtkm::cont::DataSet outData = clip.Execute(inData);

2.7.4. Density Estimation

Density estimation takes a collection of samples and estimates the density of the samples in each part of the domain (or estimate the probabilty that a sample would be at a location in the domain). The domain of samples could be a physical space, such as with particle density, or in an abstract place, such as with a histogram. The vtkm::filter::density_estimate module contains filters that estimate density in a variety of ways.

2.7.4.1. Histogram

The vtkm::filter::density_estimate::Histogram filter computes a histogram of a given scalar field.

class Histogram : public vtkm::filter::Filter

Construct the histogram of a given field.

The range of the field is evenly split to a set number of bins (set by SetNumberOfBins()). This filter then counts the number of values in the filter that are in each bin.

The result of this filter is stored in a vtkm::cont::DataSet with no points or cells. It contains only a single field containing the histogram (bin counts). The field has an association of vtkm::cont::Field::Association::WholeDataSet. The field contains an array of vtkm::Id with the bin counts. By default, the field is named “histogram”, but that can be changed with the SetOutputFieldName() method.

If this filter is run on a partitioned data set, the result will be a vtkm::cont::PartitionedDataSet containing a single vtkm::cont::DataSet as previously described.

Public Functions

inline void SetNumberOfBins(vtkm::Id count)

Set the number of bins for the resulting histogram.

By default, a histogram with 10 bins is created.

inline vtkm::Id GetNumberOfBins() const

Get the number of bins for the resulting histogram.

inline void SetRange(const vtkm::Range &range)

Set the range to use to generate the histogram.

If range is set to empty, the field’s global range (computed using vtkm::cont::FieldRangeGlobalCompute) will be used.

inline const vtkm::Range &GetRange() const

Get the range used to generate the histogram.

If the returned range is empty, then the field’s global range will be used.

inline vtkm::Float64 GetBinDelta() const

Returns the size of bin in the computed histogram.

This value is only valid after a call to Execute.

inline vtkm::Range GetComputedRange() const

Returns the range used for most recent execute.

If SetRange is used to specify a non-empty range, then this range will be returned. Otherwise, the coputed range is returned. This value is only valid after a call to Execute.

2.7.4.2. Particle Density

VTK‑m provides multiple filters to take as input a collection of points and build a regular mesh containing an estimate of the density of particles in that space. These filters inhert from vtkm::filter::density_estimate::ParticleDensityBase.

class ParticleDensityBase : public vtkm::filter::Filter

Subclassed by vtkm::filter::density_estimate::ParticleDensityCloudInCell, vtkm::filter::density_estimate::ParticleDensityNearestGridPoint

Public Functions

inline void SetComputeNumberDensity(bool flag)

Toggles between summing mass and computing instances.

When this flag is false (the default), the active field of the input is accumulated in each bin of the output. When this flag is set to true, the active field is ignored and the associated particles are simply counted.

inline bool GetComputeNumberDensity() const

Toggles between summing mass and computing instances.

When this flag is false (the default), the active field of the input is accumulated in each bin of the output. When this flag is set to true, the active field is ignored and the associated particles are simply counted.

inline void SetDivideByVolume(bool flag)

Specifies whether the accumulated mass (or count) is divided by the volume of the cell.

When this flag is on (the default), the computed mass will be divided by the volume of the bin to give a density value. Turning off this flag provides an accumulated mass or count.

inline bool GetDivideByVolume() const

Specifies whether the accumulated mass (or count) is divided by the volume of the cell.

When this flag is on (the default), the computed mass will be divided by the volume of the bin to give a density value. Turning off this flag provides an accumulated mass or count.

inline void SetDimension(const vtkm::Id3 &dimension)

The number of bins in the grid used as regions to estimate density.

To estimate particle density, this filter defines a uniform grid in space.

The numbers specify the number of bins (i.e. cells in the output mesh) in each dimension, not the number of points in the output mesh.

inline vtkm::Id3 GetDimension() const

The number of bins in the grid used as regions to estimate density.

To estimate particle density, this filter defines a uniform grid in space.

The numbers specify the number of bins (i.e. cells in the output mesh) in each dimension, not the number of points in the output mesh.

inline void SetOrigin(const vtkm::Vec3f &origin)

The lower-left (minimum) corner of the domain of density estimation.

inline vtkm::Vec3f GetOrigin() const

The lower-left (minimum) corner of the domain of density estimation.

inline void SetSpacing(const vtkm::Vec3f &spacing)

The spacing of the grid points used to form the grid for density estimation.

inline vtkm::Vec3f GetSpacing() const

The spacing of the grid points used to form the grid for density estimation.

inline void SetBounds(const vtkm::Bounds &bounds)

The bounds of the region where density estimation occurs.

This method can be used in place of SetOrigin and SetSpacing. It is often easiest to compute the bounds of the input coordinate system (or other spatial region) to use as the input.

The dimensions must be set before the bounds are set. Calling SetDimension will change the ranges of the bounds.

2.7.4.2.1. Nearest Grid Point

The vtkm::filter::density_estimate::ParticleDensityNearestGridPoint filter defines a 3D grid of bins. It then takes from the input a collection of particles, identifies which bin each particle lies in, and sums some attribute from a field of the input (or the particles can simply be counted).

class ParticleDensityNearestGridPoint : public vtkm::filter::density_estimate::ParticleDensityBase

Estimate the density of particles using the Nearest Grid Point method.

This filter takes a collection of particles. The particles are infinitesimal in size with finite mass (or other scalar attributes such as charge). The filter estimates density by imposing a regular grid (as specified by SetDimensions, SetOrigin, and SetSpacing) and summing the mass of particles within each cell in the grid. Each input particle is assigned to one bin that it falls in.

The mass of particles is established by setting the active field (using SetActiveField). Note that the “mass” can actually be another quantity. For example, you could use electrical charge in place of mass to compute the charge density. Once the sum of the mass is computed for each grid cell, the mass is divided by the volume of the cell. Thus, the density will be computed as the units of the mass field per the cubic units of the coordinate system. If you just want a sum of the mass in each cell, turn off the DivideByVolume feature of this filter. In addition, you can also simply count the number of particles in each cell by calling SetComputeNumberDensity(true).

This operation is helpful in the analysis of particle-based simulation where the data often requires conversion or deposition of particles’ attributes, such as mass, to an overlaying mesh. This allows further identification of regions of interest based on the spatial distribution of particles attributes, for example, high density regions could be considered as clusters or halos while low density regions could be considered as bubbles or cavities in the particle data.

Since there is no specific vtkm::cont::CellSet for particles in VTK-m, this filter treats the vtkm::cont::CoordinateSystem of the vtkm::cont::DataSet as the positions of the particles while ignoring the details of the vtkm::cont::CellSet.

2.7.4.2.2. Cloud in Cell

The vtkm::filter::density_estimate::ParticleDensityCloudInCell filter defines a 3D grid of bins. It then takes from the input a collection of particles, identifies which bin each particle lies in, and then redistributes each particle’s attribute to the 8 vertices of the containing bin. The filter then sums up all the contributions of particles for each bin in the grid.

class ParticleDensityCloudInCell : public vtkm::filter::density_estimate::ParticleDensityBase

Estimate the density of particles using the Cloud-in-Cell method.

This filter takes a collection of particles. The particles are infinitesimal in size with finite mass (or other scalar attributes such as charge). The filter estimates density by imposing a regular grid (as specified by SetDimensions, SetOrigin, and SetSpacing) and summing the mass of particles within each cell in the grid. The particle’s mass is divided among the 8 nearest neighboring bins. This differs from ParticleDensityNearestGridPoint, which just finds the nearest containing bin.

The mass of particles is established by setting the active field (using SetActiveField). Note that the “mass” can actually be another quantity. For example, you could use electrical charge in place of mass to compute the charge density. Once the sum of the mass is computed for each grid cell, the mass is divided by the volume of the cell. Thus, the density will be computed as the units of the mass field per the cubic units of the coordinate system. If you just want a sum of the mass in each cell, turn off the DivideByVolume feature of this filter. In addition, you can also simply count the number of particles in each cell by calling SetComputeNumberDensity(true).

This operation is helpful in the analysis of particle-based simulation where the data often requires conversion or deposition of particles’ attributes, such as mass, to an overlaying mesh. This allows further identification of regions of interest based on the spatial distribution of particles attributes, for example, high density regions could be considered as clusters or halos while low density regions could be considered as bubbles or cavities in the particle data.

2.7.4.3. Statistics

Simple descriptive statics for data in field arrays can be computed with vtkm::filter::density_estimate::Statistics.

class Statistics : public vtkm::filter::Filter

Computes descriptive statistics of an input field.

This filter computes the following statistics on the active field of the input.

  • N

  • Min

  • Max

  • Sum

  • Mean

  • M2

  • M3

  • M4

  • SampleStddev

  • PopulationStddev

  • SampleVariance

  • PopulationVariance

  • Skewness

  • Kurtosis

M2, M3, and M4 are the second, third, and fourth moments, respectively.

Note that this filter treats the “sample” and the “population” as the same with the same mean. The difference between the two forms of variance is how they are normalized. The population variance is normalized by dividing the second moment by N. The sample variance uses Bessel’s correction and divides the second moment by N-1 instead. The standard deviation, which is just the square root of the variance, follows the same difference.

The result of this filter is stored in a vtkm::cont::DataSet with no points or cells. It contains only fields with the same names as the list above. All fields have an association of vtkm::cont::Field::Association::WholeDataSet.

If Execute is called with a vtkm::cont::PartitionedDataSet, then the partitions of the output will match those of the input. Additionally, the containing vtkm::cont::PartitionedDataSet will contain the same fields associated with vtkm::cont::Field::Association::Global that provide the overall statistics of all partitions.

If this filter is used inside of an MPI job, then each vtkm::cont::DataSet result will be local to the MPI rank. If Execute is called with a vtkm::cont::PartitionedDataSet, then the fields attached to the vtkm::cont::PartitionedDataSet container will have the overall statistics across all MPI ranks (in addition to all partitions). Global MPI statistics for a single vtkm::cont::DataSet can be computed by creating a vtkm::cont::PartitionedDataSet with that as a single partition.

2.7.5. Entity Extraction

VTK‑m contains a collection of filters that extract a portion of one vtkm::cont::DataSet and construct a new vtkm::cont::DataSet based on that portion of the geometry. These filters are collected in the vtkm::filter::entity_extraction module.

2.7.5.1. External Faces

vtkm::filter::entity_extraction::ExternalFaces is a filter that extracts all the external faces from a polyhedral data set. An external face is any face that is on the boundary of a mesh. Thus, if there is a hole in a volume, the boundary of that hole will be considered external. More formally, an external face is one that belongs to only one cell in a mesh.

class ExternalFaces : public vtkm::filter::Filter

Extract external faces of a geometry.

ExternalFaces is a filter that extracts all external faces from a data set. An external face is defined is defined as a face/side of a cell that belongs only to one cell in the entire mesh.

Public Functions

inline virtual bool CanThread() const override

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.

inline bool GetCompactPoints() const

Option to remove unused points and compact result int a smaller array.

When CompactPoints is on, instead of copying the points and point fields from the input, the filter will create new compact fields without the unused elements. When off (the default), unused points will remain listed in the topology, but point fields and coordinate systems will be shallow-copied to the output.

inline void SetCompactPoints(bool value)

Option to remove unused points and compact result int a smaller array.

When CompactPoints is on, instead of copying the points and point fields from the input, the filter will create new compact fields without the unused elements. When off (the default), unused points will remain listed in the topology, but point fields and coordinate systems will be shallow-copied to the output.

inline bool GetPassPolyData() const

Specify how polygonal data (polygons, lines, and vertices) will be handled.

When on (the default), these cells will be passed to the output. When off, these cells will be removed from the output. (Because they have less than 3 topological dimensions, they are not considered to have any “faces.”)

void SetPassPolyData(bool value)

Specify how polygonal data (polygons, lines, and vertices) will be handled.

When on (the default), these cells will be passed to the output. When off, these cells will be removed from the output. (Because they have less than 3 topological dimensions, they are not considered to have any “faces.”)

2.7.5.2. Extract Geometry

The vtkm::filter::entity_extraction::ExtractGeometry filter extracts all of the cells in a vtkm::cont::DataSet that is inside or outside of an implicit function. Implicit functions are described in Chapter 2.12 (Implicit Functions). They define a function in 3D space that follow a geometric shape. The inside of the implicit function is the region of negative values.

class ExtractGeometry : public vtkm::filter::Filter

Extract a subset of geometry based on an implicit function.

Extracts from its input geometry all cells that are either completely inside or outside of a specified implicit function. Any type of data can be input to this filter.

To use this filter you must specify an implicit function. You must also specify whether to extract cells laying inside or outside of the implicit function. (The inside of an implicit function is the negative values region.) An option exists to extract cells that are neither inside or outside (i.e., boundary).

This differs from vtkm::filter::contour::ClipWithImplicitFunction in that vtkm::filter::contour::ClipWithImplicitFunction will subdivide boundary cells into new cells whereas this filter will not, producing a more “crinkly” output.

Public Functions

inline void SetImplicitFunction(const vtkm::ImplicitFunctionGeneral &func)

Specifies the implicit function to be used to perform extract geometry.

Only a limited number of implicit functions are supported. See vtkm::ImplicitFunctionGeneral for information on which ones.

inline bool GetExtractInside() const

Specify the region of the implicit function to keep cells.

Determines whether to extract the geometry that is on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline void SetExtractInside(bool value)

Specify the region of the implicit function to keep cells.

Determines whether to extract the geometry that is on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline void ExtractInsideOn()

Specify the region of the implicit function to keep cells.

Determines whether to extract the geometry that is on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline void ExtractInsideOff()

Specify the region of the implicit function to keep cells.

Determines whether to extract the geometry that is on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline bool GetExtractBoundaryCells() const

Specify whether cells on the boundary should be extracted.

The implicit function used to extract geometry is likely to intersect some of the cells of the input. If this flag is true, then any cells intersected by the implicit function are extracted and included in the output. This flag is false by default.

inline void SetExtractBoundaryCells(bool value)

Specify whether cells on the boundary should be extracted.

The implicit function used to extract geometry is likely to intersect some of the cells of the input. If this flag is true, then any cells intersected by the implicit function are extracted and included in the output. This flag is false by default.

inline void ExtractBoundaryCellsOn()

Specify whether cells on the boundary should be extracted.

The implicit function used to extract geometry is likely to intersect some of the cells of the input. If this flag is true, then any cells intersected by the implicit function are extracted and included in the output. This flag is false by default.

inline void ExtractBoundaryCellsOff()

Specify whether cells on the boundary should be extracted.

The implicit function used to extract geometry is likely to intersect some of the cells of the input. If this flag is true, then any cells intersected by the implicit function are extracted and included in the output. This flag is false by default.

inline bool GetExtractOnlyBoundaryCells() const

Specify whether to extract cells only on the boundary.

When this flag is off (the default), this filter extract the geometry in the region specified by the implicit function. When this flag is on, then only those cells that intersect the surface of the implicit function are extracted.

inline void SetExtractOnlyBoundaryCells(bool value)

GetExtractOnlyBoundaryCells.

inline void ExtractOnlyBoundaryCellsOn()

GetExtractOnlyBoundaryCells.

inline void ExtractOnlyBoundaryCellsOff()

GetExtractOnlyBoundaryCells.

2.7.5.3. Extract Points

The vtkm::filter::entity_extraction::ExtractPoints filter behaves the same as vtkm::filter::entity_extraction::ExtractGeometry (Section 2.7.5.2) except that the geometry is converted into a point cloud. The filter determines whether each point is inside or outside the implicit function and passes only those that match the criteria. The cell information of the input is thrown away and replaced with a cell set of “vertex” cells, one per point.

class ExtractPoints : public vtkm::filter::Filter

Extract only points from a geometry using an implicit function.

Extract only the points that are either inside or outside of a VTK-m implicit function. Examples include planes, spheres, boxes, etc.

Note that while any geometry type can be provided as input, the output is represented by an explicit representation of points using vtkm::cont::CellSetSingleType with one vertex cell per point.

Public Functions

inline bool GetCompactPoints() const

Option to remove unused points and compact result int a smaller array.

When CompactPoints is on, instead of copying the points and point fields from the input, the filter will create new compact fields without the unused elements. When off (the default), unused points will remain listed in the topology, but point fields and coordinate systems will be shallow-copied to the output.

inline void SetCompactPoints(bool value)

Option to remove unused points and compact result int a smaller array.

When CompactPoints is on, instead of copying the points and point fields from the input, the filter will create new compact fields without the unused elements. When off (the default), unused points will remain listed in the topology, but point fields and coordinate systems will be shallow-copied to the output.

inline void SetImplicitFunction(const vtkm::ImplicitFunctionGeneral &func)

Specifies the implicit function to be used to perform extract points.

Only a limited number of implicit functions are supported. See vtkm::ImplicitFunctionGeneral for information on which ones.

inline bool GetExtractInside() const

Specify the region of the implicit function to keep points.

Determines whether to extract the points that are on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline void SetExtractInside(bool value)

Specify the region of the implicit function to keep points.

Determines whether to extract the points that are on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline void ExtractInsideOn()

Specify the region of the implicit function to keep points.

Determines whether to extract the points that are on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

inline void ExtractInsideOff()

Specify the region of the implicit function to keep points.

Determines whether to extract the points that are on the inside of the implicit function (where the function is less than 0) or the outside (where the function is greater than 0). This flag is true by default (i.e., the interior of the implicit function will be extracted).

2.7.5.4. Extract Structured

vtkm::filter::entity_extraction::ExtractStructured is a filter that extracts a volume of interest (VOI) from a structured data set. In addition the filter is able to subsample the VOI while doing the extraction. The input and output of this filter are a structured data sets.

class ExtractStructured : public vtkm::filter::Filter

Select a piece (e.g., volume of interest) and/or subsample structured points dataset.

Select or subsample a portion of an input structured dataset. The selected portion of interested is referred to as the Volume Of Interest, or VOI. The output of this filter is a structured dataset. The filter treats input data of any topological dimension (i.e., point, line, plane, or volume) and can generate output data of any topological dimension.

To use this filter set the VOI ivar which are i-j-k min/max indices that specify a rectangular region in the data. (Note that these are 0-offset.) You can also specify a sampling rate to subsample the data.

Typical applications of this filter are to extract a slice from a volume for image processing, subsampling large volumes to reduce data size, or extracting regions of a volume with interesting data.

Public Functions

inline vtkm::RangeId3 GetVOI() const

Specifies what volume of interest (VOI) should be extracted by the filter.

The VOI is specified using the 3D indices of the structured mesh. Meshes with fewer than 3 dimensions will ignore the extra dimensions in the VOI. The VOI is inclusive on the minium index and exclusive on the maximum index.

By default the VOI is the entire input.

inline void SetVOI(vtkm::Id i0, vtkm::Id i1, vtkm::Id j0, vtkm::Id j1, vtkm::Id k0, vtkm::Id k1)

Specifies what volume of interest (VOI) should be extracted by the filter.

The VOI is specified using the 3D indices of the structured mesh. Meshes with fewer than 3 dimensions will ignore the extra dimensions in the VOI. The VOI is inclusive on the minium index and exclusive on the maximum index.

By default the VOI is the entire input.

inline void SetVOI(vtkm::Id extents[6])

Specifies what volume of interest (VOI) should be extracted by the filter.

The VOI is specified using the 3D indices of the structured mesh. Meshes with fewer than 3 dimensions will ignore the extra dimensions in the VOI. The VOI is inclusive on the minium index and exclusive on the maximum index.

By default the VOI is the entire input.

inline void SetVOI(vtkm::Id3 minPoint, vtkm::Id3 maxPoint)

Specifies what volume of interest (VOI) should be extracted by the filter.

The VOI is specified using the 3D indices of the structured mesh. Meshes with fewer than 3 dimensions will ignore the extra dimensions in the VOI. The VOI is inclusive on the minium index and exclusive on the maximum index.

By default the VOI is the entire input.

inline void SetVOI(const vtkm::RangeId3 &voi)

Specifies what volume of interest (VOI) should be extracted by the filter.

The VOI is specified using the 3D indices of the structured mesh. Meshes with fewer than 3 dimensions will ignore the extra dimensions in the VOI. The VOI is inclusive on the minium index and exclusive on the maximum index.

By default the VOI is the entire input.

inline vtkm::Id3 GetSampleRate() const

Specifies the sample rate of the VOI.

The input data can be subsampled by selecting every n-th value. The sampling can be different in each dimension. The default sampling rate is (1,1,1), meaning that no subsampling will occur.

inline void SetSampleRate(vtkm::Id i, vtkm::Id j, vtkm::Id k)

Specifies the sample rate of the VOI.

The input data can be subsampled by selecting every n-th value. The sampling can be different in each dimension. The default sampling rate is (1,1,1), meaning that no subsampling will occur.

inline void SetSampleRate(vtkm::Id3 sampleRate)

Specifies the sample rate of the VOI.

The input data can be subsampled by selecting every n-th value. The sampling can be different in each dimension. The default sampling rate is (1,1,1), meaning that no subsampling will occur.

2.7.5.5. Ghost Cell Removal

The vtkm::filter::entity_extraction::GhostCellRemove filter is used to remove cells from a data set according to a cell centered field that specifies whether a cell is a regular cell or a ghost cell. By default, the filter will get the ghost cell information that is registered in the input vtkm::cont::DataSet, but it also possible to specify an arbitrary field for this purpose.

class GhostCellRemove : public vtkm::filter::Filter

Removes cells marked as ghost cells.

This filter inspects the ghost cell field of the input and removes any cells marked as ghost cells. Although this filter nominally operates on ghost cells, other classifications, such as blanked cells, can also be recorded in the ghost cell array. See vtkm::CellClassification for the list of flags typical in a ghost array.

By default, if the input is a structured data set the filter will attempt to output a structured data set. This will be the case if all the cells along a boundary are marked as ghost cells together, which is common. If creating a structured data set is not possible, an explicit data set is produced.

Public Functions

inline void SetRemoveGhostField(bool flag)

Specify whether the ghost cell array should be removed from the input.

If this flag is true, then the ghost cell array will not be passed to the output.

inline bool GetRemoveGhostField() const

Specify whether the ghost cell array should be removed from the input.

If this flag is true, then the ghost cell array will not be passed to the output.

inline void SetTypesToRemove(vtkm::UInt8 typeFlags)

Specify which types of cells to remove.

The types to remove are specified by the flags in vtkm::CellClassification. Any cell with a ghost array flag matching one or more of these flags will be removed.

inline vtkm::UInt8 GetTypesToRemove() const

Specify which types of cells to remove.

The types to remove are specified by the flags in vtkm::CellClassification. Any cell with a ghost array flag matching one or more of these flags will be removed.

inline void SetTypesToRemoveToAll()

Set filter to remove any special cell type.

This method sets the state to remove any cell that does not have a “normal” ghost cell value of 0. Any other value represents a cell that is placeholder or otherwise not really considered part of the cell set.

inline bool AreAllTypesRemoved() const

Returns true if all abnormal cell types are removed.

inline bool GetUseGhostCellsAsField() const

Specify whether the marked ghost cells or a named field should be used as the ghost field.

When this flag is true (the default), the filter will get from the input vtkm::cont::DataSet the field (with the GetGhostCellField method). When this flag is false, the SetActiveField method of this class should be used to select which field to use as ghost cells.

inline void SetUseGhostCellsAsField(bool flag)

Specify whether the marked ghost cells or a named field should be used as the ghost field.

When this flag is true (the default), the filter will get from the input vtkm::cont::DataSet the field (with the GetGhostCellField method). When this flag is false, the SetActiveField method of this class should be used to select which field to use as ghost cells.

2.7.5.6. Threshold

A threshold operation removes topology elements from a data set that do not meet a specified criterion. The vtkm::filter::entity_extraction::Threshold filter removes all cells where the a field is outside a range of values.

Note that vtkm::filter::entity_extraction::Threshold either passes an entire cell or discards an entire cell. This can consequently lead to jagged surfaces at the interface of the threshold caused by the shape of cells that jut inside or outside the removed region. See Section 2.7.3.3 (Clip with Field) for a clipping filter that will clip off a smooth region of the mesh.

class Threshold : public vtkm::filter::Filter

Extracts cells that satisfy a threshold criterion.

Extracts all cells from any dataset type that satisfy a threshold criterion. The output of this filter stores its connectivity in a vtkm::cont::CellSetExplicit<> regardless of the input dataset type or which cells are passed.

You can threshold either on point or cell fields. If thresholding on point fields, you must specify whether a cell should be kept if some but not all of its incident points meet the criteria.

Although Threshold is primarily designed for scalar fields, there is support for thresholding on 1 or all of the components in a vector field. See the SetComponentToTest(), SetComponentToTestToAny(), and SetComponentToTestToAll() methods for more information.

Use SetActiveField() and related methods to set the field to threshold on.

Public Functions

inline void SetLowerThreshold(vtkm::Float64 value)

Specifies the lower scalar value.

Any cells where the scalar field is less than this value are removed.

inline void SetUpperThreshold(vtkm::Float64 value)

Specifies the upper scalar value.

Any cells where the scalar field is more than this value are removed.

inline vtkm::Float64 GetLowerThreshold() const

Specifies the lower scalar value.

Any cells where the scalar field is less than this value are removed.

inline vtkm::Float64 GetUpperThreshold() const

Specifies the upper scalar value.

Any cells where the scalar field is more than this value are removed.

void SetThresholdBelow(vtkm::Float64 value)

Sets the threshold criterion to pass any value less than or equal to value.

void SetThresholdAbove(vtkm::Float64 value)

Sets the threshold criterion to pass any value greater than or equal to value.

void SetThresholdBetween(vtkm::Float64 value1, vtkm::Float64 value2)

Set the threshold criterion to pass any value between (inclusive) the given values.

This method is equivalent to calling SetLowerThreshold(value1) and SetUpperThreshold(value2).

inline void SetComponentToTest(vtkm::IdComponent component)

Specifies that the threshold criteria should be applied to a specific vector component.

When thresholding on a vector field (which has more than one component per entry), the Threshold filter will by default compare the threshold criterion to the first component of the vector (component index 0). Use this method to change the component to test against.

inline void SetComponentToTestToAny()

Specifies that the threshold criteria should be applied to a specific vector component.

This method sets that the threshold criteria should be applied to all the components of the input vector field and a cell will pass if any the components match.

inline void SetComponentToTestToAll()

Specifies that the threshold criteria should be applied to a specific vector component.

This method sets that the threshold criteria should be applied to all the components of the input vector field and a cell will pass if all the components match.

inline void SetAllInRange(bool value)

Specify criteria for cells that have some points matching.

When thresholding on a point field, each cell must consider the multiple values associated with all incident points. When this flag is false (the default), the cell is passed if any of the incident points matches the threshold criterion. When this flag is true, the cell is passed only if all the incident points match the threshold criterion.

inline bool GetAllInRange() const

Specify criteria for cells that have some points matching.

When thresholding on a point field, each cell must consider the multiple values associated with all incident points. When this flag is false (the default), the cell is passed if any of the incident points matches the threshold criterion. When this flag is true, the cell is passed only if all the incident points match the threshold criterion.

inline void SetInvert(bool value)

Inverts the threshold result.

When set to true, the threshold result is inverted. That is, cells that would have been in the output with this option set to false (the default) are excluded while cells that would have been excluded from the output are included.

inline bool GetInvert() const

Inverts the threshold result.

When set to true, the threshold result is inverted. That is, cells that would have been in the output with this option set to false (the default) are excluded while cells that would have been excluded from the output are included.

2.7.6. Field Conversion

Field conversion modifies a field of a vtkm::cont::DataSet to have roughly equivalent values but with a different structure. These filters allow the field to be used in places where they otherwise would not be applicable.

2.7.6.1. Cell Average

vtkm::filter::field_conversion::CellAverage is the cell average filter. It will take a data set with a collection of cells and a field defined on the points of the data set and create a new field defined on the cells. The values of this new derived field are computed by averaging the values of the input field at all the incident points. This is a simple way to convert a point field to a cell field.

class CellAverage : public vtkm::filter::Filter

Point to cell interpolation filter.

CellAverage is a filter that transforms point data (i.e., data specified at cell points) into cell data (i.e., data specified per cell). The method of transformation is based on averaging the data values of all points used by particular cell.

The point field to convert comes from the active scalars. The default name for the output cell field is the same name as the input point field. The name can be overridden as always using the SetOutputFieldName() method.

2.7.6.2. Point Average

vtkm::filter::field_conversion::PointAverage is the point average filter. It will take a data set with a collection of cells and a field defined on the cells of the data set and create a new field defined on the points. The values of this new derived field are computed by averaging the values of the input field at all the incident cells. This is a simple way to convert a cell field to a point field.

class PointAverage : public vtkm::filter::Filter

Cell to Point interpolation filter.

PointAverage is a filter that transforms cell data (i.e., data specified per cell) into point data (i.e., data specified at cell points). The method of transformation is based on averaging the data values of all cells using a particular point.

The cell field to convert comes from the active scalars. The default name for the output cell field is the same name as the input point field. The name can be overridden as always using the SetOutputFieldName() method.

2.7.7. Field Transform

VTK‑m provides multiple filters to convert fields through some mathematical relationship.

2.7.7.1. Composite Vectors

The vtkm::filter::field_transform::CompositeVectors filter allows you to group multiple scalar fields into a single vector field. This is convenient when importing data from a souce that stores vector components in separate arrays.

class CompositeVectors : public vtkm::filter::Filter

Combine multiple scalar fields into a single vector field.

Scalar fields are selected as the active input fields, and the combined vector field is set at the output. The SetFieldNameList() method takes a std::vector of field names to use as the component fields. Alternately, the SetActiveField() method can be used to select the fields independently.

All of the input fields must be scalar values. The type of the first field determines the type of the output vector field.

Public Functions

void SetFieldNameList(const std::vector<std::string> &fieldNameList, vtkm::cont::Field::Association association = vtkm::cont::Field::Association::Any)

Specifies the names of the fields to use as components for the output.

vtkm::IdComponent GetNumberOfFields() const

The number of fields specified as inputs.

This will be the number of components in the generated field.

2.7.7.2. Cylindrical Coordinate System Transform

The vtkm::filter::field_transform::CylindricalCoordinateTransform filter is a coordinate system transformation. The filter will take a data set and transform the points of the coordinate system. By default, the filter will transform the coordinates from a Cartesian coordinate system to a cylindrical coordinate system. The order for cylindrical coordinates is \((R, \theta, Z)\). The output coordinate system will be set to the new computed coordinates.

class CylindricalCoordinateTransform : public vtkm::filter::Filter

Transform coordinates between Cartesian and cylindrical.

By default, this filter will transform the first coordinate system, but this can be changed by setting the active field.

The resulting transformation will be set as the first coordinate system in the output.

Public Functions

inline void SetCartesianToCylindrical()

Establish a transformation from Cartesian to cylindrical coordinates.

inline void SetCylindricalToCartesian()

Establish a transformation from cylindrical to Cartesian coordiantes.

2.7.7.3. Field to Colors

The vtkm::filter::field_transform::FieldToColors filter takes a field in a data set, looks up each value in a color table, and writes the resulting colors to a new field. The color to be used for each field value is specified using a vtkm::cont::ColorTable object. vtkm::cont::ColorTable objects are also used with VTK‑m’s rendering module and are described in Section 2.8.8 (Color Tables).

vtkm::filter::field_transform::FieldToColors has three modes it can use to select how it should treat the input field. These input modes are contained in vtkm::filter::field_transform::FieldToColors::InputMode. Additionally, vtkm::filter::field_transform::FieldToColors has different modes in which it can represent colors in its output. These output modes are contained in vtkm::filter::field_transform::FieldToColors::OutputMode.

class FieldToColors : public vtkm::filter::Filter

Convert an arbitrary field to an RGB or RGBA field.

This filter is useful for generating colors that could be used for rendering or other purposes.

Public Types

enum class InputMode

Identifiers used to specify how FieldToColors should treat its input scalars.

Values:

enumerator Scalar

Treat the field as a scalar field.

It is an error to provide a field of any type that cannot be directly converted to a basic floating point number (such as a vector).

enumerator Magnitude

Map the magnitude of the field.

Given a vector field, the magnitude of each field value is taken before looking it up in the color table.

enumerator Component

Map a component of a vector field as if it were a scalar.

Given a vector field, a particular component is looked up in the color table as if that component were in a scalar field. The component to map is selected with SetMappingComponent().

enum class OutputMode

Identifiers used to specify what output FieldToColors will generate.

Values:

enumerator RGB

Write out RGB fixed precision color values.

Output colors are represented as RGB values with each component represented by an unsigned byte. Specifically, these are vtkm::Vec3ui_8 values.

enumerator RGBA

Write out RGBA fixed precision color values.

Output colors are represented as RGBA values with each component represented by an unsigned byte. Specifically, these are vtkm::Vec4ui_8 values.

Public Functions

inline void SetColorTable(const vtkm::cont::ColorTable &table)

Specifies the vtkm::cont::ColorTable object to use to map field values to colors.

inline const vtkm::cont::ColorTable &GetColorTable() const

Specifies the vtkm::cont::ColorTable object to use to map field values to colors.

inline void SetMappingMode(InputMode mode)

Specify the mapping mode.

inline void SetMappingToScalar()

Treat the field as a scalar field.

It is an error to provide a field of any type that cannot be directly converted to a basic floating point number (such as a vector).

inline void SetMappingToMagnitude()

Map the magnitude of the field.

Given a vector field, the magnitude of each field value is taken before looking it up in the color table.

inline void SetMappingToComponent()

Map a component of a vector field as if it were a scalar.

Given a vector field, a particular component is looked up in the color table as if that component were in a scalar field. The component to map is selected with SetMappingComponent().

inline InputMode GetMappingMode() const

Specify the mapping mode.

inline bool IsMappingScalar() const

Returns true if this filter is in scalar mapping mode.

inline bool IsMappingMagnitude() const

Returns true if this filter is in magnitude mapping mode.

inline bool IsMappingComponent() const

Returns true if this filter is vector component mapping mode.

inline void SetMappingComponent(vtkm::IdComponent comp)

Specifies the component of the vector to use in the mapping.

This only has an effect if the input mapping mode is set to FieldToColors::InputMode::Component.

inline vtkm::IdComponent GetMappingComponent() const

Specifies the component of the vector to use in the mapping.

This only has an effect if the input mapping mode is set to FieldToColors::InputMode::Component.

inline void SetOutputMode(OutputMode mode)

Specify the output mode.

inline void SetOutputToRGB()

Write out RGB fixed precision color values.

Output colors are represented as RGB values with each component represented by an unsigned byte. Specifically, these are vtkm::Vec3ui_8 values.

inline void SetOutputToRGBA()

Write out RGBA fixed precision color values.

Output colors are represented as RGBA values with each component represented by an unsigned byte. Specifically, these are vtkm::Vec4ui_8 values.

inline OutputMode GetOutputMode() const

Specify the output mode.

inline bool IsOutputRGB() const

Returns true if this filter is in RGB output mode.

inline bool IsOutputRGBA() const

Returns true if this filter is in RGBA output mode.

void SetNumberOfSamplingPoints(vtkm::Int32 count)

Specifies how many samples to use when looking up color values.

The implementation of FieldToColors first builds an array of color samples to quickly look up colors for particular values. The size of this lookup array can be adjusted with this parameter. By default, an array of 256 colors is used.

inline vtkm::Int32 GetNumberOfSamplingPoints() const

Specifies how many samples to use when looking up color values.

The implementation of FieldToColors first builds an array of color samples to quickly look up colors for particular values. The size of this lookup array can be adjusted with this parameter. By default, an array of 256 colors is used.

2.7.7.4. Generate Ids

The vtkm::filter::field_transform::GenerateIds filter creates point and/or cell fields that mimic the identifier for the respective element.

class GenerateIds : public vtkm::filter::Filter

Adds fields to a vtkm::cont::DataSet that give the ids for the points and cells.

This filter will add (by default) a point field named pointids that gives the index of the associated point and likewise a cell field named cellids for the associated cell indices. These fields are useful for tracking the provenance of the elements of a vtkm::cont::DataSet as it gets manipulated by filters. It is also convenient for adding indices to operations designed for fields and generally creating test data.

Public Functions

inline const std::string &GetPointFieldName() const

The name given to the generated point field.

By default, the name is pointids.

inline void SetPointFieldName(const std::string &name)

The name given to the generated point field.

By default, the name is pointids.

inline const std::string &GetCellFieldName() const

The name given to the generated cell field.

By default, the name is cellids.

inline void SetCellFieldName(const std::string &name)

The name given to the generated cell field.

By default, the name is cellids.

inline bool GetGeneratePointIds() const

Specify whether the point id field is generated.

When GeneratePointIds is true (the default), a field echoing the point indices is generated. When set to false, this output is not created.

inline void SetGeneratePointIds(bool flag)

Specify whether the point id field is generated.

When GeneratePointIds is true (the default), a field echoing the point indices is generated. When set to false, this output is not created.

inline bool GetGenerateCellIds() const

Specify whether the cell id field is generated.

When GenerateCellIds is true (the default), a field echoing the cell indices is generated. When set to false, this output is not created.

inline void SetGenerateCellIds(bool flag)

Specify whether the cell id field is generated.

When GenerateCellIds is true (the default), a field echoing the cell indices is generated. When set to false, this output is not created.

inline bool GetUseFloat() const

Specify whether the generated fields should be integer or float.

When UseFloat is false (the default), then the fields generated will have type vtkm::Id. If it is set to true, then the fields will be generated with type vtkm::FloatDefault.

inline void SetUseFloat(bool flag)

Specify whether the generated fields should be integer or float.

When UseFloat is false (the default), then the fields generated will have type vtkm::Id. If it is set to true, then the fields will be generated with type vtkm::FloatDefault.

2.7.7.5. Log Values

The vtkm::filter::field_transform::LogValues filter can be used to take the logarithm of all values in a field. The filter is able to take the logarithm to a number of predefined bases identified by vtkm::filter::field_transform::LogValues::LogBase.

class LogValues : public vtkm::filter::Filter

Adds field to a vtkm::cont::DataSet that gives the log values for the user specified field.

By default, LogValues takes a natural logarithm (of base e). The base of the logarithm can be set to one of the bases listed in LogBase with SetBaseValue().

Logarithms are often used to rescale data to simultaneously show data at different orders of magnitude. It allows small changes in small numbers be visible next to much larger numbers with less precision. One problem with this approach is if there exist numbers very close to zero, the scale at the low range could make all but the smallest numbers comparatively hard to see. Thus, LogValues supports setting a minimum value (with SetMinValue()) that will clamp any smaller values to that.

Public Types

enum class LogBase

Identifies a type of logarithm as specified by its base.

Values:

enumerator E

Take the natural logarithm.

The logarithm is set to the mathematical constant e (about 2.718). This is a constant that has many uses in calculus and other mathematics, and a logarithm of base e is often referred to as the “natural” logarithm.

enumerator TWO

Take the base 2 logarithm.

The base 2 logarithm is particularly useful for estimating the depth of a binary hierarchy.

enumerator TEN

Take the base 10 logarithm.

The base 10 logarithm is handy to convert a number to its order of magnitude based on our standard base 10 human counting system.

Public Functions

inline const LogBase &GetBaseValue() const

Specify the base of the logarithm.

inline void SetBaseValue(const LogBase &base)

Specify the base of the logarithm.

inline void SetBaseValueToE()

Take the natural logarithm.

The logarithm is set to the mathematical constant e (about 2.718). This is a constant that has many uses in calculus and other mathematics, and a logarithm of base e is often referred to as the “natural” logarithm.

inline void SetBaseValueTo2()

Take the base 2 logarithm.

The base 2 logarithm is particularly useful for estimating the depth of a binary hierarchy.

inline void SetBaseValueTo10()

Take the base 10 logarithm.

The base 10 logarithm is handy to convert a number to its order of magnitude based on our standard base 10 human counting system.

inline vtkm::FloatDefault GetMinValue() const

Specifies the minimum value to take the logarithm of.

Before taking the logarithm, this filter will check the value to this minimum value and clamp it to the minimum value if it is lower. This is useful to prevent values from approching negative infinity.

By default, no minimum value is used.

inline void SetMinValue(const vtkm::FloatDefault &value)

Specifies the minimum value to take the logarithm of.

Before taking the logarithm, this filter will check the value to this minimum value and clamp it to the minimum value if it is lower. This is useful to prevent values from approching negative infinity.

By default, no minimum value is used.

2.7.7.6. Point Elevation

The vtkm::filter::field_transform::PointElevation filter computes the “elevation” of a field of point coordinates in space. Example 2.19 gives a demonstration of the elevation filter.

class PointElevation : public vtkm::filter::Filter

Generate a scalar field along a specified direction.

The filter will take a data set and a field of 3 dimensional vectors and compute the distance along a line defined by a low point and a high point. Any point in the plane touching the low point and perpendicular to the line is set to the minimum range value in the elevation whereas any point in the plane touching the high point and perpendicular to the line is set to the maximum range value. All other values are interpolated linearly between these two planes. This filter is commonly used to compute the elevation of points in some direction, but can be repurposed for a variety of measures.

The default name for the output field is `elevation’, but that can be overridden as always using the SetOutputFieldName() method.

Public Functions

inline void SetLowPoint(const vtkm::Vec3f_64 &point)

Specify the coordinate of the low point.

The plane of low values is defined by the plane that contains the low point and is normal to the direction from the low point to the high point. All vector values on this plane are assigned the low value.

inline void SetLowPoint(vtkm::Float64 x, vtkm::Float64 y, vtkm::Float64 z)

Specify the coordinate of the low point.

The plane of low values is defined by the plane that contains the low point and is normal to the direction from the low point to the high point. All vector values on this plane are assigned the low value.

inline void SetHighPoint(const vtkm::Vec3f_64 &point)

Specify the coordinate of the high point.

The plane of high values is defined by the plane that contains the high point and is normal to the direction from the low point to the high point. All vector values on this plane are assigned the high value.

inline void SetHighPoint(vtkm::Float64 x, vtkm::Float64 y, vtkm::Float64 z)

Specify the coordinate of the high point.

The plane of high values is defined by the plane that contains the high point and is normal to the direction from the low point to the high point. All vector values on this plane are assigned the high value.

inline void SetRange(vtkm::Float64 low, vtkm::Float64 high)

Specify the range of values to output.

Values at the low plane are given low and values at the high plane are given high. Values in between the planes have a linearly interpolated value based on the relative distance between the two planes.

2.7.7.7. Point Transform

The vtkm::filter::field_transform::PointTransform filter performs affine transforms is the point transform filter.

class PointTransform : public vtkm::filter::Filter

Perform affine transforms to point coordinates or vector fields.

This filter will take a data set and a field of 3 dimensional vectors and perform the specified point transform operation. Several methods are provided to apply many common affine transformations (e.g., translation, rotation, and scale). You can also provide a general 4x4 transformation matrix with SetTransform().

The main use case for PointTransform is to perform transformations of objects in 3D space, which is done by applying these transforms to the coordinate system. This filter will operate on the vtkm::cont::CoordinateSystem of the input data unless a different active field is specified. Likewise, this filter will save its results as the first coordinate system in the output unless SetChangeCoordinateSystem() is set to say otherwise.

The default name for the output field is “transform”, but that can be overridden as always using the SetOutputFieldName() method.

Public Functions

inline void SetTranslation(const vtkm::FloatDefault &tx, const vtkm::FloatDefault &ty, const vtkm::FloatDefault &tz)

Translates, or moves, each point in the input field by a given direction.

inline void SetTranslation(const vtkm::Vec3f &v)

Translates, or moves, each point in the input field by a given direction.

inline void SetRotation(const vtkm::FloatDefault &angleDegrees, const vtkm::Vec3f &axis)

Rotate the input field about a given axis.

Parameters:
  • angleDegrees[in] The amount of rotation to perform, given in degrees.

  • axis[in] The rotation is made around a line that goes through the origin and pointing in this direction in the counterclockwise direction.

inline void SetRotation(const vtkm::FloatDefault &angleDegrees, const vtkm::FloatDefault &axisX, const vtkm::FloatDefault &axisY, const vtkm::FloatDefault &axisZ)

Rotate the input field about a given axis.

The rotation is made around a line that goes through the origin and pointing in the direction specified by axisX, axisY, and axisZ in the counterclockwise direction.

Parameters:
  • angleDegrees[in] The amount of rotation to perform, given in degrees.

  • axisX[in] The X value of the rotation axis.

  • axisY[in] The Y value of the rotation axis.

  • axisZ[in] The Z value of the rotation axis.

inline void SetRotationX(const vtkm::FloatDefault &angleDegrees)

Rotate the input field around the X axis by the given degrees.

inline void SetRotationY(const vtkm::FloatDefault &angleDegrees)

Rotate the input field around the Y axis by the given degrees.

inline void SetRotationZ(const vtkm::FloatDefault &angleDegrees)

Rotate the input field around the Z axis by the given degrees.

inline void SetScale(const vtkm::FloatDefault &s)

Scale the input field.

Each coordinate is multiplied by tghe associated scale factor.

inline void SetScale(const vtkm::FloatDefault &sx, const vtkm::FloatDefault &sy, const vtkm::FloatDefault &sz)

Scale the input field.

Each coordinate is multiplied by tghe associated scale factor.

inline void SetScale(const vtkm::Vec3f &v)

Scale the input field.

Each coordinate is multiplied by tghe associated scale factor.

inline void SetTransform(const vtkm::Matrix<vtkm::FloatDefault, 4, 4> &mtx)

Set a general transformation matrix.

Each field value is multiplied by this 4x4 as a homogeneous coordinate. That is a 1 component is added to the end of each 3D vector to put it in the form [x, y, z, 1]. The matrix is then premultiplied to this as a column vector.

The functions in vtkm/Transform3D.h can be used to help build these transform matrices.

void SetChangeCoordinateSystem(bool flag)

Specify whether the result should become the coordinate system of the output.

When this flag is on (the default) the first coordinate system in the output vtkm::cont::DataSet is set to the transformed point coordinates.

2.7.7.8. Spherical Coordinate System Transform

The vtkm::filter::field_transform::SphericalCoordinateTransform filter is a coordinate system transformation. The filter will take a data set and transform the points of the coordinate system. By default, the filter will transform the coordinates from a Cartesian coordinate system to a spherical coordinate system. The order for spherical coordinates is \((R, \theta, \phi)\) where \(R\) is the radius, \(\theta\) is the azimuthal angle and \(\phi\) is the polar angle. The output coordinate system will be set to the new computed coordinates.

class SphericalCoordinateTransform : public vtkm::filter::Filter

Transform coordinates between Cartesian and spherical.

By default, this filter will transform the first coordinate system, but this can be changed by setting the active field.

The resulting transformation will be set as the first coordinate system in the output.

Public Functions

inline void SetCartesianToSpherical()

Establish a transformation from Cartesian to spherical coordinates.

inline void SetSphericalToCartesian()

Establish a transformation from spherical to Cartesian coordiantes.

2.7.7.9. Warp

The vtkm::filter::field_transform::Warp filter modifies points in a vtkm::cont::DataSet by moving points along scaled direction vectors. By default, the vtkm::filter::field_transform::Warp filter modifies the coordinate system and writes its results to the coordiante system. A vector field can be selected as directions, or a constant direction can be specified. A constant direction is particularly useful for generating a carpet plot. A scalar field can be selected to scale the displacement, and a constant scale factor adjustment can be specified.

class Warp : public vtkm::filter::Filter

Modify points by moving points along scaled direction vectors.

This filter displaces the point coordinates of a dataset either in the direction of a direction vector field or in a constant direction.

The filter starts with a set of point coordinates or other vectors. By default these vectors are the coordinate system, but they can be changed by modifying active field 0. These vectors are then displaced by a set of vectors. This is done by selecting a field of directions, a field of scales, and an additional scale factor. The directions are multiplied by the scale field and the scale factor, and this displacement is added to the vector.

It is common to wish to warp in a constant direction by a scaled amount. To support this so called “WarpScalar”, the Warp filter allows you to specify a constant direction direction with the SetConstantDirection() method. When this is set, no direction field is retrieved. By default Warp uses (0, 0, 1) as the direction direction.

It is also common to wish to simply apply a vector direction field (with a possible constant scale). To support this so called “WarpVector”, the Warp filter allows you to ignore the scale field with the SetUseScaleField() method. When this is unset, no scale field is retrieved. Calling SetScaleField() turns on the UseScaleField flag. By default, Warp uses will not use the scale field unless specified.

The main use case for Warp is to adjust the spatial location and shape of objects in 3D space. This filter will operate on the vtkm::cont::CoordinateSystem of the input data unless a different active field is specified. Likewise, this filter will save its results as the first coordinate system in the output unless SetChangeCoordinateSystem() is set to say otherwise.

Subclassed by vtkm::filter::field_transform::WarpScalar, vtkm::filter::field_transform::WarpVector

Public Functions

inline void SetDirectionField(const std::string &name)

Specify a field to use as the directions.

The directions, when not set to use constant directions, are set as active field index 1.

inline std::string GetDirectionFieldName() const

Specify a field to use as the directions.

The directions, when not set to use constant directions, are set as active field index 1.

inline void SetConstantDirection(const vtkm::Vec3f &direction)

Specify a constant value to use as the directions.

This will provide a (constant) direction of the direction, and the direction field will be ignored.

inline const vtkm::Vec3f &GetConstantDirection() const

Specify a constant value to use as the directions.

This will provide a (constant) direction of the direction, and the direction field will be ignored.

inline void SetUseConstantDirection(bool flag)

Specifies whether a direction field or a constant direction direction is used.

When true, the constant direction direction is used. When false, the direction field (active field index 1) is used.

inline bool GetUseConstantDirection() const

Specifies whether a direction field or a constant direction direction is used.

When true, the constant direction direction is used. When false, the direction field (active field index 1) is used.

inline void SetScaleField(const std::string &name)

Specify a field to use to scale the directions.

The scale factor field scales the size of the direction. The scale factor, when not set to use a constant factor, is set as active field index 2.

inline std::string GetScaleFieldName() const

Specify a field to use to scale the directions.

The scale factor field scales the size of the direction. The scale factor, when not set to use a constant factor, is set as active field index 2.

inline void SetUseScaleField(bool flag)

Specifies whether a scale factor field is used.

When true, a scale factor field the constant scale factor is used. When false, the scale factor field (active field index 2) is used.

inline bool GetUseScaleField() const

Specifies whether a scale factor field is used.

When true, a scale factor field the constant scale factor is used. When false, the scale factor field (active field index 2) is used.

inline void SetScaleFactor(vtkm::FloatDefault scale)

Specifies an additional scale factor to scale the displacements.

When using a non-constant scale field, it is possible that the scale field is of the wrong units and needs to be rescaled. This scale factor is multiplied to the direction and scale to re-adjust the overall scale.

inline vtkm::FloatDefault GetScaleFactor() const

Specifies an additional scale factor to scale the displacements.

When using a non-constant scale field, it is possible that the scale field is of the wrong units and needs to be rescaled. This scale factor is multiplied to the direction and scale to re-adjust the overall scale.

inline void SetChangeCoordinateSystem(bool flag)

Specify whether the result should become the coordinate system of the output.

When this flag is on (the default) the first coordinate system in the output vtkm::cont::DataSet is set to the transformed point coordinates.

inline bool GetChangeCoordinateSystem() const

Specify whether the result should become the coordinate system of the output.

When this flag is on (the default) the first coordinate system in the output vtkm::cont::DataSet is set to the transformed point coordinates.

2.7.8. Flow Analysis

Flow visualization is used to analyze vector fields that represent the movement of a fluid. The basic operation of most flow visualization algorithms is particle advection, which traces the path a particle would take given the direction and speed dictated by the vector field. There are multiple ways in which to represent flow in this manner, and consequently VTK‑m contains several filters that trace streams in different ways. These filters inherit from vtkm::filter::flow::FilterParticleAdvection, which provides several important methods.

class FilterParticleAdvection : public vtkm::filter::Filter

base class for advecting particles in a vector field.

Takes as input a vector field and seed locations and advects the seeds through the flow field.

Subclassed by vtkm::filter::flow::FilterParticleAdvectionSteadyState< ParticleAdvection >, vtkm::filter::flow::FilterParticleAdvectionSteadyState< WarpXStreamline >, vtkm::filter::flow::FilterParticleAdvectionSteadyState< Streamline >, vtkm::filter::flow::FilterParticleAdvectionUnsteadyState< PathParticle >, vtkm::filter::flow::FilterParticleAdvectionUnsteadyState< Pathline >, vtkm::filter::flow::FilterParticleAdvectionSteadyState< Derived >, vtkm::filter::flow::FilterParticleAdvectionUnsteadyState< Derived >

Public Functions

inline virtual bool CanThread() const override

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.

inline void SetStepSize(vtkm::FloatDefault s)

Specifies the step size used for the numerical integrator.

The numerical integrators operate by advancing each particle by a finite amount. This parameter defines the distance to advance each time. Smaller values are more accurate but take longer to integrate. An appropriate step size is usually around the size of each cell.

inline void SetNumberOfSteps(vtkm::Id n)

Specifies the maximum number of integration steps for each particle.

Some particle paths may loop and continue indefinitely. This parameter sets an upper limit on the total length of advection.

template<typename ParticleType>
inline void SetSeeds(vtkm::cont::ArrayHandle<ParticleType> &seeds)

Specify the seed locations for the particle advection.

Each seed represents one particle that is advected by the vector field. The particles are represented by a vtkm::Particle object or similar type of object (such as vtkm::ChargedParticle).

template<typename ParticleType>
inline void SetSeeds(const std::vector<ParticleType> &seeds, vtkm::CopyFlag copyFlag = vtkm::CopyFlag::On)

Specify the seed locations for the particle advection.

Each seed represents one particle that is advected by the vector field. The particles are represented by a vtkm::Particle object or similar type of object (such as vtkm::ChargedParticle).

Flow filters operate either on a “steady state” flow that does not change or on an “unsteady state” flow that is continually changing over time. An unsteady state filter must be executed multiple times for subsequent time steps. The filter operates with data from two adjacent time steps. This is managed by the vtkm::filter::flow::FilterParticleAdvectionUnsteadyState superclass.

2.7.8.1. Streamlines

Streamlines are a powerful technique for the visualization of flow fields. A streamline is a curve that is parallel to the velocity vector of the flow field. Individual streamlines are computed from an initial point location (seed) using a numerical method to integrate the point through the flow field.

class Streamline : public vtkm::filter::flow::FilterParticleAdvectionSteadyState<Streamline>

Advect particles in a vector field and display the path they take.

This filter takes as input a velocity vector field and seed locations. It then traces the path each seed point would take if moving at the velocity specified by the field. Mathematically, this is the curve that is tangent to the velocity field everywhere.

The output of this filter is a vtkm::cont::DataSet containing a collection of poly-lines representing the paths the seed particles take.

The vtkm::filter::flow::Streamline filter also uses several inherited methods: vtkm::filter::flow::FilterParticleAdvection::SetSeeds(), vtkm::filter::flow::FilterParticleAdvection::SetStepSize(), and vtkm::filter::flow::FilterParticleAdvection::SetNumberOfSteps().

 1  vtkm::filter::flow::Streamline streamlines;
 2
 3  // Specify the seeds.
 4  vtkm::cont::ArrayHandle<vtkm::Particle> seedArray;
 5  seedArray.Allocate(2);
 6  seedArray.WritePortal().Set(0, vtkm::Particle({ 0, 0, 0 }, 0));
 7  seedArray.WritePortal().Set(1, vtkm::Particle({ 1, 1, 1 }, 1));
 8
 9  streamlines.SetActiveField("vectorvar");
10  streamlines.SetStepSize(0.1f);
11  streamlines.SetNumberOfSteps(100);
12  streamlines.SetSeeds(seedArray);
13
14  vtkm::cont::DataSet output = streamlines.Execute(inData);

2.7.8.2. Pathlines

Pathlines are the analog to streamlines for time varying vector fields. Individual pathlines are computed from an initial point location (seed) using a numerical method to integrate the point through the flow field.

This filter requires two data sets as input, which represent the data for two sequential time steps. The “Previous” data set, which marks the data at the earlier time step, is passed into the filter throught the standard Execute method. The “Next” data set, which marks the data at the later time step, is specified as state to the filter using methods.

class Pathline : public vtkm::filter::flow::FilterParticleAdvectionUnsteadyState<Pathline>

Advect particles in a time-varying vector field and display the path they take.

This filter takes as input a velocity vector field, changing between two time steps, and seed locations. It then traces the path each seed point would take if moving at the velocity specified by the field.

The output of this filter is a vtkm::cont::DataSet containing a collection of poly-lines representing the paths the seed particles take.

As an unsteady state flow filter, vtkm::filter::flow::Pathline must be executed multiple times for subsequent time steps. The filter operates with data from two adjacent time steps. This is managed by the vtkm::filter::flow::FilterParticleAdvectionUnsteadyState superclass.

The vtkm::filter::flow::Pathline filter uses several other inherited methods: vtkm::filter::flow::FilterParticleAdvectionUnsteadyState::SetPreviousTime(), vtkm::filter::flow::FilterParticleAdvectionUnsteadyState::SetNextTime(), vtkm::filter::flow::FilterParticleAdvectionUnsteadyState::SetNextDataSet(), vtkm::filter::flow::FilterParticleAdvection::SetSeeds(), vtkm::filter::flow::FilterParticleAdvection::SetStepSize(), and vtkm::filter::flow::FilterParticleAdvection::SetNumberOfSteps().

 1  vtkm::filter::flow::Pathline pathlines;
 2
 3  // Specify the seeds.
 4  vtkm::cont::ArrayHandle<vtkm::Particle> seedArray;
 5  seedArray.Allocate(2);
 6  seedArray.WritePortal().Set(0, vtkm::Particle({ 0, 0, 0 }, 0));
 7  seedArray.WritePortal().Set(1, vtkm::Particle({ 1, 1, 1 }, 1));
 8
 9  pathlines.SetActiveField("vectorvar");
10  pathlines.SetStepSize(0.1f);
11  pathlines.SetNumberOfSteps(100);
12  pathlines.SetSeeds(seedArray);
13  pathlines.SetPreviousTime(0.0f);
14  pathlines.SetNextTime(1.0f);
15  pathlines.SetNextDataSet(inData2);
16
17  vtkm::cont::DataSet pathlineCurves = pathlines.Execute(inData1);

2.7.8.3. Stream Surface

A stream surface is defined as a continuous surface that is everywhere tangent to a specified vector field. The vtkm::filter::flow::StreamSurface filter computes a stream surface from a set of input points and the vector field of the input data set. The stream surface is created by creating streamlines from each input point and then connecting adjacent streamlines with a series of triangles.

class StreamSurface : public vtkm::filter::Filter

Generate stream surfaces from a vector field.

This filter takes as input a velocity vector field and seed locations. The seed locations should be arranged in a line or curve. The filter then traces the path each seed point would take if moving at the velocity specified by the field and connects all the lines together into a surface. Mathematically, this is the surface that is tangent to the velocity field everywhere.

The output of this filter is a vtkm::cont::DataSet containing a mesh for the created surface.

Public Functions

inline void SetStepSize(vtkm::FloatDefault s)

Specifies the step size used for the numerical integrator.

The numerical integrators operate by advancing each particle by a finite amount. This parameter defines the distance to advance each time. Smaller values are more accurate but take longer to integrate. An appropriate step size is usually around the size of each cell.

inline void SetNumberOfSteps(vtkm::Id n)

Specifies the maximum number of integration steps for each particle.

Some particle paths may loop and continue indefinitely. This parameter sets an upper limit on the total length of advection.

template<typename ParticleType>
inline void SetSeeds(vtkm::cont::ArrayHandle<ParticleType> &seeds)

Specify the seed locations for the particle advection.

Each seed represents one particle that is advected by the vector field. The particles are represented by a vtkm::Particle object.

template<typename ParticleType>
inline void SetSeeds(const std::vector<ParticleType> &seeds, vtkm::CopyFlag copyFlag = vtkm::CopyFlag::On)

Specify the seed locations for the particle advection.

Each seed represents one particle that is advected by the vector field. The particles are represented by a vtkm::Particle object.

 1  vtkm::filter::flow::StreamSurface streamSurface;
 2
 3  // Specify the seeds.
 4  vtkm::cont::ArrayHandle<vtkm::Particle> seedArray;
 5  seedArray.Allocate(2);
 6  seedArray.WritePortal().Set(0, vtkm::Particle({ 0, 0, 0 }, 0));
 7  seedArray.WritePortal().Set(1, vtkm::Particle({ 1, 1, 1 }, 1));
 8
 9  streamSurface.SetActiveField("vectorvar");
10  streamSurface.SetStepSize(0.1f);
11  streamSurface.SetNumberOfSteps(100);
12  streamSurface.SetSeeds(seedArray);
13
14  vtkm::cont::DataSet output = streamSurface.Execute(inData);

2.7.8.4. Lagrangian Coherent Structures

Lagrangian coherent structures (LCS) are distinct structures present in a flow field that have a major influence over nearby trajectories over some interval of time. Some of these structures may be sources, sinks, saddles, or vortices in the flow field. Identifying Lagrangian coherent structures is part of advanced flow analysis and is an important part of studying flow fields. These structures can be studied by calculating the finite time Lyapunov exponent (FTLE) for a flow field at various locations, usually over a regular grid encompassing the entire flow field. If the provided input dataset is structured, then by default the points in this data set will be used as seeds for advection. The vtkm::filter::flow::LagrangianStructures filter is used to compute the FTLE of a flow field.

class LagrangianStructures : public vtkm::filter::Filter

Compute the finite time Lyapunov exponent (FTLE) of a vector field.

The FTLE is computed by advecting particles throughout the vector field and analyizing where they diverge or converge. By default, the points of the input vtkm::cont::DataSet are all advected for this computation unless an auxiliary grid is established.

Public Functions

inline virtual bool CanThread() const override

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.

inline void SetStepSize(vtkm::FloatDefault s)

Specifies the step size used for the numerical integrator.

The numerical integrators operate by advancing each particle by a finite amount. This parameter defines the distance to advance each time. Smaller values are more accurate but take longer to integrate. An appropriate step size is usually around the size of each cell.

inline vtkm::FloatDefault GetStepSize()

Specifies the step size used for the numerical integrator.

The numerical integrators operate by advancing each particle by a finite amount. This parameter defines the distance to advance each time. Smaller values are more accurate but take longer to integrate. An appropriate step size is usually around the size of each cell.

inline void SetNumberOfSteps(vtkm::Id n)

Specify the maximum number of steps each particle is allowed to traverse.

This can limit the total length of displacements used when computing the FTLE.

inline vtkm::Id GetNumberOfSteps()

Specify the maximum number of steps each particle is allowed to traverse.

This can limit the total length of displacements used when computing the FTLE.

inline void SetAdvectionTime(vtkm::FloatDefault advectionTime)

Specify the time interval for the advection.

The FTLE works by advecting all points a finite distance, and this parameter specifies how far to advect.

inline vtkm::FloatDefault GetAdvectionTime()

Specify the time interval for the advection.

The FTLE works by advecting all points a finite distance, and this parameter specifies how far to advect.

inline void SetUseAuxiliaryGrid(bool useAuxiliaryGrid)

Specify whether to use an auxiliary grid.

When this flag is off (the default), then the points of the mesh representing the vector field are advected and used for computing the FTLE. However, if the mesh is too coarse, the FTLE will likely be inaccurate. Or if the mesh is unstructured the FTLE may be less efficient to compute. When this flag is on, an auxiliary grid of uniformly spaced points is used for the FTLE computation.

inline bool GetUseAuxiliaryGrid()

Specify whether to use an auxiliary grid.

When this flag is off (the default), then the points of the mesh representing the vector field are advected and used for computing the FTLE. However, if the mesh is too coarse, the FTLE will likely be inaccurate. Or if the mesh is unstructured the FTLE may be less efficient to compute. When this flag is on, an auxiliary grid of uniformly spaced points is used for the FTLE computation.

inline void SetAuxiliaryGridDimensions(vtkm::Id3 auxiliaryDims)

Specify the dimensions of the auxiliary grid for FTLE calculation.

Seeds for advection will be placed along the points of this auxiliary grid. This option has no effect unless the UseAuxiliaryGrid option is on.

inline vtkm::Id3 GetAuxiliaryGridDimensions()

Specify the dimensions of the auxiliary grid for FTLE calculation.

Seeds for advection will be placed along the points of this auxiliary grid. This option has no effect unless the UseAuxiliaryGrid option is on.

inline void SetUseFlowMapOutput(bool useFlowMapOutput)

Specify whether to use flow maps instead of advection.

If the start and end points for FTLE calculation are known already, advection is an unnecessary step. This flag allows users to bypass advection, and instead use a precalculated flow map. By default this option is off.

inline bool GetUseFlowMapOutput()

Specify whether to use flow maps instead of advection.

If the start and end points for FTLE calculation are known already, advection is an unnecessary step. This flag allows users to bypass advection, and instead use a precalculated flow map. By default this option is off.

inline void SetOutputFieldName(std::string outputFieldName)

Specify the name of the output field in the data set returned.

By default, the field will be named FTLE.

inline std::string GetOutputFieldName()

Specify the name of the output field in the data set returned.

By default, the field will be named FTLE.

inline void SetFlowMapOutput(vtkm::cont::ArrayHandle<vtkm::Vec3f> &flowMap)

Specify the array representing the flow map output to be used for FTLE calculation.

inline vtkm::cont::ArrayHandle<vtkm::Vec3f> GetFlowMapOutput()

Specify the array representing the flow map output to be used for FTLE calculation.

2.7.9. Geometry Refinement

Geometry refinement modifies the geometry of a vtkm::cont::DataSet. It might add, change, or remove components of the structure, but the general representation will be the same.

2.7.9.1. Convert to a Point Cloud

Data in a vtkm::cont::DataSet is typically connected together by cells in a mesh structure. However, it is sometimes the case where data are simply represented as a cloud of unconnected points. These meshless data sets are best represented in a vtkm::cont::DataSet by a collection of “vertex” cells.

The vtkm::filter::geometry_refinement::ConvertToPointCloud filter converts a data to a point cloud. It does this by throwing away any existing cell set and replacing it with a collection of vertex cells, one per point. vtkm::filter::geometry_refinement::ConvertToPointCloud is useful to add a cell set to a vtkm::cont::DataSet that has points but no cells. It is also useful to treat data as a collection of sample points rather than an interconnected mesh.

class ConvertToPointCloud : public vtkm::filter::Filter

Convert a DataSet to a point cloud.

A point cloud in VTK-m is represented as a data set with “vertex” shape cells. This filter replaces the CellSet in a DataSet with a CellSet of only vertex cells. There will be one cell per point.

This filter is useful for dropping the cells of any DataSet so that you can operate on it as just a collection of points. It is also handy for completing a DataSet that does not have a CellSet associated with it or has points that do not belong to cells.

Note that all fields associated with cells are dropped. This is because the cells are dropped.

Public Functions

inline void SetAssociateFieldsWithCells(bool flag)

By default, all the input point fields are kept as point fields in the output.

However, the output has exactly one cell per point and it might be easier to treat the fields as cell fields. When this flag is turned on, the point field association is changed to cell.

Note that any field that is marked as point coordinates will remain as point fields. It is not valid to set a cell field as the point coordinates.

inline bool GetAssociateFieldsWithCells() const

By default, all the input point fields are kept as point fields in the output.

However, the output has exactly one cell per point and it might be easier to treat the fields as cell fields. When this flag is turned on, the point field association is changed to cell.

Note that any field that is marked as point coordinates will remain as point fields. It is not valid to set a cell field as the point coordinates.

2.7.9.2. Shrink

The vtkm::filter::geometry_refinement::Shrink independently reduces the size of each class. Rather than uniformly reduce the size of the whole data set (which can be done with vtkm::filter::field_transform::PointTransform), this filter separates the cells from each other and shrinks them around their centroid. This is useful for making an “exploded view” of the data where the facets of the data are moved away from each other to see inside.

class Shrink : public vtkm::filter::Filter

Shrink cells of an arbitrary dataset by a constant factor.

The Shrink filter shrinks the cells of a DataSet towards their centroid, computed as the average position of the cell points. This filter disconnects the cells, duplicating the points connected to multiple cells. The resulting CellSet is always an ExplicitCellSet.

Public Functions

inline void SetShrinkFactor(vtkm::FloatDefault factor)

Specify the scale factor to size each cell.

The shrink factor specifies the ratio of the shrunk cell to its original size. This value must be between 0 and 1. A value of 1 is the same size as the input, and a value of 0 shrinks each cell to a point.

inline vtkm::FloatDefault GetShrinkFactor() const

Specify the scale factor to size each cell.

The shrink factor specifies the ratio of the shrunk cell to its original size. This value must be between 0 and 1. A value of 1 is the same size as the input, and a value of 0 shrinks each cell to a point.

2.7.9.3. Split Sharp Edges

The vtkm::filter::geometry_refinement::SplitSharpEdges filter splits sharp manifold edges where the feature angle between the adjacent surfaces are larger than a threshold value. This is most useful to preserve sharp edges when otherwise applying smooth shading during rendering.

class SplitSharpEdges : public vtkm::filter::Filter

Split sharp polygon mesh edges with a large feature angle between the adjacent cells.

Split sharp manifold edges where the feature angle between the adjacent polygonal cells are larger than a threshold value. The feature angle is the angle between the normals of the two polygons. Two polygons on the same plane have a feature angle of 0. Perpendicular polygons have a feature angle of 90 degrees.

When an edge is split, it adds a new point to the coordinates and updates the connectivity of an adjacent surface. For example, consider two adjacent triangles (0,1,2) and (2,1,3) where edge (1,2) needs to be split. Two new points 4 (duplication of point 1) and 5 (duplication of point 2) would be added and the later triangle’s connectivity would be changed to (5,4,3). By default, all old point’s fields would be copied to the new point.

Note that “split” edges do not have space added between them. They are still adjacent visually, but the topology becomes disconnectered there. Splitting sharp edges is most useful to duplicate normal shading vectors to make a sharp shading effect.

Public Functions

inline void SetFeatureAngle(vtkm::FloatDefault value)

Specify the feature angle threshold to split on.

The feature angle is the angle between the normals of the two polygons. Two polygons on the same plane have a feature angle of 0. Perpendicular polygons have a feature angle of 90 degrees.

Any edge with a feature angle larger than this threshold will be split. The feature angle is specified in degrees. The default value is 30 degrees.

inline vtkm::FloatDefault GetFeatureAngle() const

Specify the feature angle threshold to split on.

The feature angle is the angle between the normals of the two polygons. Two polygons on the same plane have a feature angle of 0. Perpendicular polygons have a feature angle of 90 degrees.

Any edge with a feature angle larger than this threshold will be split. The feature angle is specified in degrees. The default value is 30 degrees.

2.7.9.4. Tetrahedralize

The vtkm::filter::geometry_refinement::Tetrahedralize filter converts all the polyhedra in a vtkm::cont::DataSet into tetrahedra.

class Tetrahedralize : public vtkm::filter::Filter

Convert all polyhedra of a vtkm::cont::DataSet into tetrahedra.

Note that although the tetrahedra will occupy the same space of the cells that they replace, the interpolation of point fields within these cells might differ. For example, the first order interpolation of a hexahedron uses trilinear interpolation, which actually results in cubic equations. This differs from the purely linear field in a tetrahedron, so the tetraheda replacement of the hexahedron will not have exactly the same interpolation.

2.7.9.5. Triangulate

The vtkm::filter::geometry_refinement::Triangulate filter converts all the polyhedra in a vtkm::cont::DataSet into tetrahedra.

class Triangulate : public vtkm::filter::Filter

Convert all polygons of a vtkm::cont::DataSet into triangles.

Note that although the triangles will occupy the same space of the cells that they replace, the interpolation of point fields within these cells might differ. For example, the first order interpolation of a quadrilateral uses bilinear interpolation, which actually results in quadratic equations. This differs from the purely linear field in a triangle, so the triangle replacement of the quadrilateral will not have exactly the same interpolation.

2.7.9.6. Tube

The vtkm::filter::geometry_refinement::Tube filter generates a tube around each line and polyline in the input data set.

class Tube : public vtkm::filter::Filter

Generate a tube around each line and polyline.

The radius, number of sides, and end capping can be specified for each tube. The orientation of the geometry of the tube are computed automatically using a heuristic to minimize the twisting along the input data set.

Public Functions

inline void SetRadius(vtkm::FloatDefault r)

Specify the radius of each tube.

inline void SetNumberOfSides(vtkm::Id n)

Specify the number of sides for each tube.

The tubes are generated using a polygonal approximation. This option determines how many facets will be generated around the tube.

inline void SetCapping(bool v)

The Tube filter can optionally add a cap at the ends of each tube.

This option specifies whether that cap is generated.

1  vtkm::filter::geometry_refinement::Tube tubeFilter;
2
3  tubeFilter.SetRadius(0.5f);
4  tubeFilter.SetNumberOfSides(7);
5  tubeFilter.SetCapping(true);
6
7  vtkm::cont::DataSet output = tubeFilter.Execute(inData);

2.7.9.7. Vertex Clustering

The vtkm::filter::geometry_refinement::VertexClustering filter simplifies a polygonal mesh. It does so by dividing space into a uniform grid of bin and then merges together all points located in the same bin. The smaller the dimensions of this binning grid, the fewer polygons will be in the output cells and the coarser the representation. This surface simplification is an important operation to support level of detail (LOD) rendering in visualization applications.

class VertexClustering : public vtkm::filter::Filter

Reduce the number of triangles in a mesh.

VertexClustering is a filter to reduce the number of triangles in a triangle mesh, forming a good approximation to the original geometry. The input must be a vtkm::cont::DataSet that contains only triangles.

The general approach of the algorithm is to cluster vertices in a uniform binning of space, accumulating to an average point within each bin. In more detail, the algorithm first gets the bounds of the input poly data. It then breaks this bounding volume into a user-specified number of spatial bins. It then reads each triangle from the input and hashes its vertices into these bins. Then, if 2 or more vertices of the triangle fall in the same bin, the triangle is dicarded. If the triangle is not discarded, it adds the triangle to the list of output triangles as a list of vertex identifiers. (There is one vertex id per bin.) After all the triangles have been read, the representative vertex for each bin is computed. This determines the spatial location of the vertices of each of the triangles in the output.

To use this filter, specify the divisions defining the spatial subdivision in the x, y, and z directions. Compared to algorithms such as vtkQuadricClustering, a significantly higher bin count is recommended as it doesn’t increase the computation or memory of the algorithm and will produce significantly better results.

Public Functions

inline void SetNumberOfDivisions(const vtkm::Id3 &num)

Specifies the dimensions of the uniform grid that establishes the bins used for clustering.

Setting smaller numbers of dimensions produces a smaller output, but with a coarser representation of the surface.

inline const vtkm::Id3 &GetNumberOfDivisions() const

Specifies the dimensions of the uniform grid that establishes the bins used for clustering.

Setting smaller numbers of dimensions produces a smaller output, but with a coarser representation of the surface.

1  vtkm::filter::geometry_refinement::VertexClustering vertexClustering;
2
3  vertexClustering.SetNumberOfDivisions(vtkm::Id3(128, 128, 128));
4
5  vtkm::cont::DataSet simplifiedSurface = vertexClustering.Execute(originalSurface);

2.7.10. Mesh Information

VTK‑m provides several filters that derive information about the structure of the geometry. This can be information about the shape of cells or their connections.

2.7.10.1. Cell Size Measurements

The vtkm::filter::mesh_info::CellMeasures filter integrates the size of each cell in a mesh and reports the size in a new cell field.

class CellMeasures : public vtkm::filter::Filter

Compute the size measure of each cell in a dataset.

CellMeasures is a filter that generates a new cell data array (i.e., one value specified per cell) holding the signed measure of the cell or 0 (if measure is not well defined or the cell type is unsupported).

By default, the new cell-data array is named “measure”.

Public Functions

inline void SetMeasure(vtkm::filter::mesh_info::IntegrationType measure)

Specify the type of integrations to support.

This filter can support integrating the size of 1D elements (arclength measurements), 2D elements (area measurements), and 3D elements (volume measurements). The measures to perform are specified with a vtkm::filter::mesh_info::IntegrationType.

By default, the size measure for all types of elements is performed.

inline vtkm::filter::mesh_info::IntegrationType GetMeasure() const

Specify the type of integrations to support.

This filter can support integrating the size of 1D elements (arclength measurements), 2D elements (area measurements), and 3D elements (volume measurements). The measures to perform are specified with a vtkm::filter::mesh_info::IntegrationType.

By default, the size measure for all types of elements is performed.

inline void SetMeasureToArcLength()

Compute the length of 1D elements.

inline void SetMeasureToArea()

Compute the area of 2D elements.

inline void SetMeasureToVolume()

Compute the volume of 3D elements.

inline void SetMeasureToAll()

Compute the size of all types of elements.

inline void SetCellMeasureName(const std::string &name)

Specify the name of the field generated.

If not set, measure is used.

inline const std::string &GetCellMeasureName() const

Specify the name of the field generated.

If not set, measure is used.

By default, vtkm::filter::mesh_info::CellMeasures will compute the measures of all types of cells. It is sometimes desirable to limit the types of cells to measure to prevent the resulting field from mixing values of different units. The appropriate measure to compute can be specified with the vtkm::filter::mesh_info::IntegrationType enumeration.

enum class vtkm::filter::mesh_info::IntegrationType

Specifies over what types of mesh elements CellMeasures will operate.

The values of IntegrationType may be |-ed together to select multiple

Values:

enumerator None
enumerator ArcLength

Compute the length of 1D elements.

enumerator Area

Compute the area of 2D elements.

enumerator Volume

Compute the volume of 3D elements.

enumerator AllMeasures

Compute the size of all types of elements.

2.7.10.2. Ghost Cell Classification

The vtkm::filter::mesh_info::GhostCellClassify filter determines which cells should be considered ghost cells in a structured data set. The ghost cells are expected to be on the border.

class GhostCellClassify : public vtkm::filter::Filter

Determines which cells should be considered ghost cells in a structured data set.

The ghost cells are expected to be on the border. The outer layer of cells are marked as ghost cells and the remainder marked as normal.

This filter generates a new cell-centered field marking the status of each cell. Each entry is set to either vtkm::CellClassification::Normal or vtkm::CellClassification::Ghost.

Public Functions

inline void SetGhostCellName(const std::string &fieldName)

Set the name of the output field name.

The output field is also marked as the ghost cell field in the output vtkm::cont::DataSet.

inline const std::string &GetGhostCellName()

Set the name of the output field name.

The output field is also marked as the ghost cell field in the output vtkm::cont::DataSet.

2.7.10.3. Mesh Quality Metrics

VTK‑m provides several filters to compute metrics about the mesh quality. These filters produce a new cell field that holds a given metric for the shape of the cell. The metrics for this filter come from the Verdict library, and full mathematical descriptions for each metric can be found in the Verdict documentation (Sandia technical report SAND2007-1751, https://coreform.com/papers/verdict_quality_library.pdf).

class MeshQualityArea : public vtkm::filter::Filter

Compute the area of each cell.

This only produces values for triangles and quadrilaterals.

Public Functions

vtkm::Float64 ComputeTotalArea(const vtkm::cont::DataSet &input)

Computes the area of all polygonal cells and returns the total area.

vtkm::Float64 ComputeAverageArea(const vtkm::cont::DataSet &input)

Computes the average area of cells.

This method first computes the total area of all cells and then divides that by the number of cells in the dataset.

class MeshQualityAspectGamma : public vtkm::filter::Filter

For each cell, compute the normalized root-mean-square of the edge lengths.

This only produces values for tetrahedra.

The root-mean-square edge length is normalized to the volume such that the value is 1 for an equilateral tetrahedron. The acceptable range for good quality meshes is considered to be [1, 3]. The normal range of values is [1, FLOAT_MAX].

class MeshQualityAspectRatio : public vtkm::filter::Filter

Compute for each cell the ratio of its longest edge to its circumradius.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

An acceptable range of this mesh for a good quality polygon is [1, 1.3], and the acceptable range for a good quality polyhedron is [1, 3]. Normal values for any cell type have the range [1, FLOAT_MAX].

class MeshQualityCondition : public vtkm::filter::Filter

Compute for each cell the condition number of the weighted Jacobian matrix.

This only produces values for triangles, quadrilaterals, and tetrahedra.

The acceptable range of values for a good quality cell is [1, 1.3] for triangles, [1, 4] for quadrilaterals, and [1, 3] for tetrahedra.

class MeshQualityDiagonalRatio : public vtkm::filter::Filter

Compute for each cell the ratio of the maximum diagonal to the minimum diagonal.

This only produces values for quadrilaterals and hexahedra.

An acceptable range for a good quality cell is [0.65, 1]. The normal range is [0, 1], but a degenerate cell with no size will have the value of infinity.

class MeshQualityDimension : public vtkm::filter::Filter

Compute for each cell a metric specifically designed for Sandia’s Pronto code.

This only produces values for hexahedra.

class MeshQualityJacobian : public vtkm::filter::Filter

Compute for each cell the minimum determinant of the Jacobian matrix, over corners and cell center.

This only produces values for quadrilaterals, tetrahedra, and hexahedra.

class MeshQualityMaxAngle : public vtkm::filter::Filter

Computes the maximum angle within each cell in degrees.

This only produces values for triangles and quadrilaterals.

For a good quality triangle, this value should be in the range [60, 90]. Poorer quality triangles can have a value as high as 180. For a good quality quadrilateral, this value should be in the range [90, 135]. Poorer quality quadrilaterals can have a value as high as 360.

class MeshQualityMaxDiagonal : public vtkm::filter::Filter

Computes the maximum diagonal length within each cell in degrees.

This only produces values for hexahedra.

class MeshQualityMinAngle : public vtkm::filter::Filter

Computes the minimum angle within each cell in degrees.

This only produces values for triangles and quadrilaterals.

For a good quality triangle, this value should be in the range [30, 60]. Poorer quality triangles can have a value as low as 0. For a good quality quadrilateral, this value should be in the range [45, 90]. Poorer quality quadrilaterals can have a value as low as 0.

class MeshQualityMinDiagonal : public vtkm::filter::Filter

Computes the minimal diagonal length within each cell in degrees.

This only produces values for hexahedra.

class MeshQualityOddy : public vtkm::filter::Filter

Compute for each cell the maximum deviation of a metric tensor from an identity matrix, over all corners and cell center.

This only produces values for quadrilaterals and hexahedra.

For a good quality quadrilateral or hexahedron, this value should be in the range [0, 0.5]. Poorer quality cells can have unboundedly larger values.

class MeshQualityRelativeSizeSquared : public vtkm::filter::Filter

Compute for each cell the ratio of area or volume to the mesh average.

If S is the size of a cell and avgS is the average cell size in the mesh, then let R = S/avgS. R is “normalized” to be in the range [0, 1] by taking the minimum of R and 1/R. This value is then squared.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For a good quality triangle, the relative sized squared should be in the range [0.25, 1]. For a good quality quadrilateral, it should be in the range [0.3, 1]. For a good quality tetrahedron, it should be in the range [0.3, 1]. For a good quality hexahedron, it should be in the range [0.5, 1]. Poorer quality cells can have a relative size squared as low as 0.

class MeshQualityScaledJacobian : public vtkm::filter::Filter

Compute for each cell a metric derived from the Jacobian matric with normalization involving edge length.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For a triangle, an acceptable range for good quality is [0.5, 2*sqrt(3)/3]. The value for an equalateral triangle is 1. The normal range is [-2*sqrt(3)/3), 2*sqrt(3)/3], but malformed cells can have plus or minus the maximum float value.

For a quadrilateral, an acceptable range for good quality is [0.3, 1]. The unit square has a value of 1. The normal range as well as the full range is [-1, 1].

For a tetrahedron, an acceptable range for good quality is [0.5, sqrt(2)/2]. The value for a unit equalateral triangle is 1. The normal range of values is [-sqrt(2)/2, sqrt(2)/2], but malformed cells can have plus or minus the maximum float value.

For a hexahedron, an acceptable range for good quality is [0.5, 1]. The unit cube has a value of 1. The normal range is [ -1, 1 ], but malformed cells can have a maximum float value.

class MeshQualityShape : public vtkm::filter::Filter

Compute a shape-based metric for each cell.

This metric is based on the condition number of the Jacobian matrix.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For good quality triangles, the acceptable range is [0.25, 1]. Good quality quadrilaterals, tetrahedra, hexahedra are in the range [0.3, 1]. Poorer quality cells can have values as low as 0.

class MeshQualityShapeAndSize : public vtkm::filter::Filter

Compute a metric for each cell based on the shape scaled by the cell size.

This filter multiplies the values of the shape metric by the relative size squared metric. See vtkm::filter::mesh_info::MeshQualityShape and vtkm::filter::mesh_info::MeshQualityRelativeSizeSquared for details on each of those metrics.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For a good quality cell, this value will be in the range [0.2, 1]. Poorer quality cells can have values as low as 0.

class MeshQualityShear : public vtkm::filter::Filter

Compute the shear of each cell.

The shear of a cell is computed by taking the minimum of the Jacobian at each corner divided by the length of the corner’s adjacent edges.

This only produces values for quadrilaterals and hexahedra. Good quality cells will have values in the range [0.3, 1]. Poorer quality cells can have values as low as 0.

class MeshQualitySkew : public vtkm::filter::Filter

Compute the skew of each cell.

The skew is computed as the dot product between unit vectors in the principal directions. (For 3D objects, the skew is taken as the maximum of all planes.)

This only produces values for quadrilaterals and hexahedra.

Good quality cells will have a skew in the range [0, 0.5]. A unit square or cube will have a skew of 0. Poor quality cells can have a skew up to 1 although a malformed cell might have its skew be infinite.

class MeshQualityStretch : public vtkm::filter::Filter

Compute the stretch of each cell.

The stretch of a cell is computed as the ratio of the minimum edge length to the maximum diagonal, normalized for the unit cube. A good quality cell will have a stretch in the range [0.25, 1]. Poorer quality cells can have a stretch as low as 0 although a malformed cell might return a strech of infinity.

This only produces values for quadrilaterals and hexahedra.

class MeshQualityTaper : public vtkm::filter::Filter

Compute the taper of each cell.

The taper of a quadrilateral is computed as the maximum ratio of the cross-derivative with its shortest associated principal axis.

This only produces values for quadrilaterals and hexahedra.

A good quality quadrilateral will have a taper in the range of [0, 0.7]. A good quality hexahedron will have a taper in the range of [0, 0.5]. The unit square or cube will have a taper of 0. Poorer quality cells will have larger values (with no upper limit).

class MeshQualityVolume : public vtkm::filter::Filter

Compute the volume each cell.

This only produces values for tetrahedra, pyramids, wedges, and hexahedra.

Public Functions

vtkm::Float64 ComputeTotalVolume(const vtkm::cont::DataSet &input)

Computes the volume of all polyhedral cells and returns the total area.

vtkm::Float64 ComputeAverageVolume(const vtkm::cont::DataSet &input)

Computes the average volume of cells.

This method first computes the total volume of all cells and then divides that by the number of cells in the dataset.

class MeshQualityWarpage : public vtkm::filter::Filter

Compute the flatness of cells.

This only produces values for quadrilaterals. It is defined as the cosine of the minimum dihedral angle formed by the planes intersecting in diagonals (to the fourth power).

This metric will be 1 for a perfectly flat quadrilateral and be lower as the quadrilateral deviates from the plane. A good quality quadrilateral will have a value in the range [0.3, 1]. Poorer quality cells having lower values down to -1, although malformed cells might have an infinite value.

Note that the value of this filter is consistent with the equivalent metric in VisIt, and it differs from the implementation in the Verdict library. The Verdict library returns 1 - value.

The vtkm::filter::mesh_info::MeshQuality filter consolidates all of these metrics into a single filter. The metric to compute is selected with the vtkm::filter::mesh_info::MeshQuality::SetMetric().

class MeshQuality : public vtkm::filter::Filter

Computes the quality of an unstructured cell-based mesh.

The quality is defined in terms of the summary statistics (frequency, mean, variance, min, max) of metrics computed over the mesh cells. One of several different metrics can be specified for a given cell type, and the mesh can consist of one or more different cell types. The resulting mesh quality is stored as one or more new fields in the output dataset of this filter, with a separate field for each cell type. Each field contains the metric summary statistics for the cell type. Summary statists with all 0 values imply that the specified metric does not support the cell type.

Public Functions

void SetMetric(CellMetric metric)

Specify the metric to compute on the mesh.

inline CellMetric GetMetric() const

Specify the metric to compute on the mesh.

std::string GetMetricName() const

Return a string describing the metric selected.

The metric to compute is identified using the vtkm::filter::mesh_info::CellMetric enum.

enum class vtkm::filter::mesh_info::CellMetric

Values:

enumerator Area

Compute the area of each cell.

This only produces values for triangles and quadrilaterals.

enumerator AspectGamma

For each cell, compute the normalized root-mean-square of the edge lengths.

This only produces values for tetrahedra.

The root-mean-square edge length is normalized to the volume such that the value is 1 for an equilateral tetrahedron. The acceptable range for good quality meshes is considered to be [1, 3]. The normal range of values is [1, FLOAT_MAX].

enumerator AspectRatio

Compute for each cell the ratio of its longest edge to its circumradius.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

An acceptable range of this mesh for a good quality polygon is [1, 1.3], and the acceptable range for a good quality polyhedron is [1, 3]. Normal values for any cell type have the range [1, FLOAT_MAX].

enumerator Condition

Compute for each cell the condition number of the weighted Jacobian matrix.

This only produces values for triangles, quadrilaterals, and tetrahedra.

The acceptable range of values for a good quality cell is [1, 1.3] for triangles, [1, 4] for quadrilaterals, and [1, 3] for tetrahedra.

enumerator DiagonalRatio

Compute for each cell the ratio of the maximum diagonal to the minimum diagonal.

This only produces values for quadrilaterals and hexahedra.

An acceptable range for a good quality cell is [0.65, 1]. The normal range is [0, 1], but a degenerate cell with no size will have the value of infinity.

enumerator Dimension

Compute for each cell a metric specifically designed for Sandia’s Pronto code.

This only produces values for hexahedra.

enumerator Jacobian

Compute for each cell the minimum determinant of the Jacobian matrix, over corners and cell center.

This only produces values for quadrilaterals, tetrahedra, and hexahedra.

enumerator MaxAngle

Computes the maximum angle within each cell in degrees.

This only produces values for triangles and quadrilaterals.

For a good quality triangle, this value should be in the range [60, 90]. Poorer quality triangles can have a value as high as 180. For a good quality quadrilateral, this value should be in the range [90, 135]. Poorer quality quadrilaterals can have a value as high as 360.

enumerator MaxDiagonal

Computes the maximum diagonal length within each cell in degrees.

This only produces values for hexahedra.

enumerator MinAngle

Computes the minimum angle within each cell in degrees.

This only produces values for triangles and quadrilaterals.

For a good quality triangle, this value should be in the range [30, 60]. Poorer quality triangles can have a value as low as 0. For a good quality quadrilateral, this value should be in the range [45, 90]. Poorer quality quadrilaterals can have a value as low as 0.

enumerator MinDiagonal

Computes the minimal diagonal length within each cell in degrees.

This only produces values for hexahedra.

enumerator Oddy

Compute for each cell the maximum deviation of a metric tensor from an identity matrix, over all corners and cell center.

This only produces values for quadrilaterals and hexahedra.

For a good quality quadrilateral or hexahedron, this value should be in the range [0, 0.5]. Poorer quality cells can have unboundedly larger values.

enumerator RelativeSizeSquared

Compute for each cell the ratio of area or volume to the mesh average.

If S is the size of a cell and avgS is the average cell size in the mesh, then let R = S/avgS. R is “normalized” to be in the range [0, 1] by taking the minimum of R and 1/R. This value is then squared.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For a good quality triangle, the relative sized squared should be in the range [0.25, 1]. For a good quality quadrilateral, it should be in the range [0.3, 1]. For a good quality tetrahedron, it should be in the range [0.3, 1]. For a good quality hexahedron, it should be in the range [0.5, 1]. Poorer quality cells can have a relative size squared as low as 0.

enumerator ScaledJacobian

Compute for each cell a metric derived from the Jacobian matric with normalization involving edge length.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For a triangle, an acceptable range for good quality is [0.5, 2*sqrt(3)/3]. The value for an equalateral triangle is 1. The normal range is [-2*sqrt(3)/3), 2*sqrt(3)/3], but malformed cells can have plus or minus the maximum float value.

For a quadrilateral, an acceptable range for good quality is [0.3, 1]. The unit square has a value of 1. The normal range as well as the full range is [-1, 1].

For a tetrahedron, an acceptable range for good quality is [0.5, sqrt(2)/2]. The value for a unit equalateral triangle is 1. The normal range of values is [-sqrt(2)/2, sqrt(2)/2], but malformed cells can have plus or minus the maximum float value.

For a hexahedron, an acceptable range for good quality is [0.5, 1]. The unit cube has a value of 1. The normal range is [ -1, 1 ], but malformed cells can have a maximum float value.

enumerator Shape

Compute a shape-based metric for each cell.

This metric is based on the condition number of the Jacobian matrix.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For good quality triangles, the acceptable range is [0.25, 1]. Good quality quadrilaterals, tetrahedra, hexahedra are in the range [0.3, 1]. Poorer quality cells can have values as low as 0.

enumerator ShapeAndSize

Compute a metric for each cell based on the shape scaled by the cell size.

This filter multiplies the values of the shape metric by the relative size squared metric. See vtkm::filter::mesh_info::MeshQualityShape and vtkm::filter::mesh_info::MeshQualityRelativeSizeSquared for details on each of those metrics.

This only produces values for triangles, quadrilaterals, tetrahedra, and hexahedra.

For a good quality cell, this value will be in the range [0.2, 1]. Poorer quality cells can have values as low as 0.

enumerator Shear

Compute the shear of each cell.

The shear of a cell is computed by taking the minimum of the Jacobian at each corner divided by the length of the corner’s adjacent edges.

This only produces values for quadrilaterals and hexahedra. Good quality cells will have values in the range [0.3, 1]. Poorer quality cells can have values as low as 0.

enumerator Skew

Compute the skew of each cell.

The skew is computed as the dot product between unit vectors in the principal directions. (For 3D objects, the skew is taken as the maximum of all planes.)

This only produces values for quadrilaterals and hexahedra.

Good quality cells will have a skew in the range [0, 0.5]. A unit square or cube will have a skew of 0. Poor quality cells can have a skew up to 1 although a malformed cell might have its skew be infinite.

enumerator Stretch

Compute the stretch of each cell.

The stretch of a cell is computed as the ratio of the minimum edge length to the maximum diagonal, normalized for the unit cube. A good quality cell will have a stretch in the range [0.25, 1]. Poorer quality cells can have a stretch as low as 0 although a malformed cell might return a strech of infinity.

This only produces values for quadrilaterals and hexahedra.

enumerator Taper

Compute the taper of each cell.

The taper of a quadrilateral is computed as the maximum ratio of the cross-derivative with its shortest associated principal axis.

This only produces values for quadrilaterals and hexahedra.

A good quality quadrilateral will have a taper in the range of [0, 0.7]. A good quality hexahedron will have a taper in the range of [0, 0.5]. The unit square or cube will have a taper of 0. Poorer quality cells will have larger values (with no upper limit).

enumerator Volume

Compute the volume each cell.

This only produces values for tetrahedra, pyramids, wedges, and hexahedra.

enumerator Warpage

Compute the flatness of cells.

This only produces values for quadrilaterals. It is defined as the cosine of the minimum dihedral angle formed by the planes intersecting in diagonals (to the fourth power).

This metric will be 1 for a perfectly flat quadrilateral and be lower as the quadrilateral deviates from the plane. A good quality quadrilateral will have a value in the range [0.3, 1]. Poorer quality cells having lower values down to -1, although malformed cells might have an infinite value.

Note that the value of this filter is consistent with the equivalent metric in VisIt, and it differs from the implementation in the Verdict library. The Verdict library returns 1 - value.

enumerator None

2.7.11. Multi-Block

Data with multiple blocks are stored in vtkm::cont::PartitionedDataSet objects. Most VTK‑m filters operate correctly on vtkm::cont::PartitionedDataSet just like they do with vtkm::cont::DataSet. However, there are some filters that are designed with operations specific to multi-block datasets.

2.7.11.1. AMR Arrays

An AMR mesh is a vtkm::cont::PartitionedDataSet with a special structure in the partitions. Each partition has a vtkm::cont::CellSetStructured cell set. The partitions form a hierarchy of grids where each level of the hierarchy refines the one above.

vtkm::cont::PartitionedDataSet does not explicitly store the structure of an AMR grid. The vtkm::filter::multi_block::AmrArrays filter determines the hierarchical structure of the AMR partitions and stores information about them in cell field arrays on each partition.

class AmrArrays : public vtkm::filter::Filter

Generate arrays describing the AMR structure in a partitioned data set.

AMR grids are represented by vtkm::cont::PartitionedDataSet, but this class does not explicitly store the hierarchical structure of the mesh refinement. This hierarchical arrangement needs to be captured in fields that describe where blocks reside in the hierarchy. This filter analyses the arrangement of partitions in a vtkm::cont::PartitionedDataSet and generates the following field arrays.

  • vtkAmrLevel The AMR level at which the partition resides (with 0 being the most coarse level). All the values for a particular partition are set to the same value.

  • vtkAmrIndex A unique identifier for each partition of a particular level. Each partition of the same level will have a unique index, but the indices will repeat across levels. All the values for a particular partition are set to the same value.

  • vtkCompositeIndex A unique identifier for each partition. This index is the same as the index used for the partition in the containing vtkm::cont::PartitionedDataSet. All the values for a particular partition are set to the same value.

  • vtkGhostType It is common for refinement levels in an AMR structure to overlap more coarse grids. In this case, the overlapped coarse cells have invalid data. The vtkGhostType field will track which cells are overlapped and should be ignored. This array will have a 0 value for all valid cells and a non-zero value for all invalid cells. (Specifically, if the bit specified by vtkm::CellClassification::BLANKED is set, then the cell is overlapped with a cell in a finer level.)

These arrays are stored as cell fields in the partitions.

This filter only operates on partitioned data sets where all the partitions have cell sets of type vtkm::cont::CellSetStructured. This is characteristic of AMR data sets.

Did You Know?

The names of the generated field arrays arrays (e.g. vtkAmrLevel) are chosen to be compatible with the equivalent arrays in VTK. This is why they use the prefix of “vtk” instead of “vtkm”. Likewise, the flags used for vtkGhostType are compatible with VTK.

2.7.11.2. Merge Data Sets

A vtkm::cont::PartitionedDataSet can often be treated the same as a vtkm::cont::DataSet as both can be passed to a filter’s Execute method. However, it is sometimes important to have all the data contained in a single DataSet. The vtkm::filter::multi_block::MergeDataSets filter can do just that to the partitions of a vtkm::cont::PartitionedDataSet.

class MergeDataSets : public vtkm::filter::Filter

Merging multiple data sets into one data set.

This filter merges multiple data sets into one data set. We assume that the input data sets have the same coordinate system. If there are missing fields in a specific data set, the filter uses the InvalidValue specified by the user to fill in the associated position of the field array.

MergeDataSets is used by passing a vtkm::cont::PartitionedDataSet to its Execute() method. The Execute() will return a vtkm::cont::PartitionedDataSet because that is the common interface for all filters. However, the vtkm::cont::PartitionedDataSet will have one partition that is all the blocks merged together.

Public Functions

inline void SetInvalidValue(vtkm::Float64 invalidValue)

Specify the value to use where field values are missing.

One issue when merging blocks in a paritioned dataset is that the blocks/partitions may have different fields. That is, one partition might not have all the fields of another partition. When these partitions are merged together, the values for this missing field must be set to something. They will be set to this value, which defaults to NaN.

inline vtkm::Float64 GetInvalidValue()

Specify the value to use where field values are missing.

One issue when merging blocks in a paritioned dataset is that the blocks/partitions may have different fields. That is, one partition might not have all the fields of another partition. When these partitions are merged together, the values for this missing field must be set to something. They will be set to this value, which defaults to NaN.

2.7.12. Resampling

All data in vtkm::cont::DataSet objects are discrete representations. It is sometimes necessary to resample this data in different ways.

2.7.12.1. Histogram Sampling

The vtkm::filter::resampling::HistSampling filter randomly samples the points of an input data set. The sampling is random but adaptive to preserve rare field value points.

class HistSampling : public vtkm::filter::Filter

Adaptively sample points to preserve tail features.

This filter randomly samples the points of a vtkm::cont::DataSet and generates a new vtkm::cont::DataSet

with a subsampling of the points. The sampling is adaptively selected to preserve tail and outlying features of the active field. That is, the more rare a field value is, the more likely the point will be selected in the sampling. This is done by creating a histogram of the field and using that to derive the importance level of each field value. Details of the algorithm can be found in the paper “In Situ Data-Driven Adaptive Sampling

for Large-scale Simulation Data Summarization” by Biswas, Dutta, Pulido, and Ahrens as published in

In Situ Infrastructures for Enabling Extreme-scale Analysis and Visualization (ISAV 2018).

The cell set of the input data is removed and replaced with a set with a vertex cell for each point. This effectively converts the data to a point cloud.

Public Functions

inline void SetNumberOfBins(vtkm::Id numberOfBins)

Specify the number of bins used when computing the histogram.

The histogram is used to select the importance of each field value. More rare field values are more likely to be selected.

inline vtkm::Id GetNumberOfBins()

Specify the number of bins used when computing the histogram.

The histogram is used to select the importance of each field value. More rare field values are more likely to be selected.

inline void SetSampleFraction(vtkm::FloatDefault fraction)

Specify the fraction of points to create in the sampled data.

A fraction of 1 means that all the points will be sampled and be in the output. A fraction of 0 means that none of the points will be sampled. A fraction of 0.5 means that half the points will be selected to be in the output.

inline vtkm::FloatDefault GetSampleFraction() const

Specify the fraction of points to create in the sampled data.

A fraction of 1 means that all the points will be sampled and be in the output. A fraction of 0 means that none of the points will be sampled. A fraction of 0.5 means that half the points will be selected to be in the output.

inline void SetSeed(vtkm::UInt32 seed)

Specify the seed used for random number generation.

The random numbers are used to select which points to pull from the input. If the same seed is used for multiple invocations, the results will be the same.

inline vtkm::UInt32 GetSeed()

Specify the seed used for random number generation.

The random numbers are used to select which points to pull from the input. If the same seed is used for multiple invocations, the results will be the same.

2.7.12.2. Probe

The vtkm::filter::resampling::Probe filter maps the fields of one vtkm::cont::DataSet onto another. This is useful for redefining meshes as well as comparing field data from two data sets with different geometries.

class Probe : public vtkm::filter::Filter

Sample the fields of a data set at specified locations.

The vtkm::filter::resampling::Probe filter samples the fields of one vtkm::cont::DataSet and places them in the fields of another vtkm::cont::DataSet.

To use this filter, first specify a geometry to probe with with SetGeometry(). The most important feature of this geometry is its coordinate system. When you call Execute(), the output will be the data specified with SetGeometry() but will have the fields of the input to Execute() transferred to it. The fields are transfered by probing the input data set at the point locations of the geometry.

Public Functions

inline void SetGeometry(const vtkm::cont::DataSet &geometry)

Specify the geometry to probe with.

When Execute() is called, the input data will be probed at all the point locations of this geometry as specified by its coordinate system.

inline const vtkm::cont::DataSet &GetGeometry() const

Specify the geometry to probe with.

When Execute() is called, the input data will be probed at all the point locations of this geometry as specified by its coordinate system.

inline void SetInvalidValue(vtkm::Float64 invalidValue)

Specify the value to use for points outside the bounds of the input.

It is possible that the sampling geometry will have points outside the bounds of the input. When this happens, the field will be set to this “invalid” value. By default, the invalid value is NaN.

inline vtkm::Float64 GetInvalidValue() const

Specify the value to use for points outside the bounds of the input.

It is possible that the sampling geometry will have points outside the bounds of the input. When this happens, the field will be set to this “invalid” value. By default, the invalid value is NaN.

2.7.13. Vector Analysis

VTK‑m’s vector analysis filters compute operations on fields related to vectors (usually in 3-space).

2.7.13.1. Cross Product

The vtkm::filter::vector_analysis::CrossProduct filter computes the cross product of two vector fields for every element in the input data set. The cross product filter computes (PrimaryField × SecondaryField). The cross product computation works for either point or cell centered vector fields.

class CrossProduct : public vtkm::filter::Filter

Compute the cross product of 3D vector fields.

The left part of the operand is the “primary” field and the right part of the operand is the “secondary” field.

Public Functions

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

Specify the primary field to operate on.

In the cross product operation A x B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline const std::string &GetPrimaryFieldName() const

Specify the primary field to operate on.

In the cross product operation A x B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::cont::Field::Association GetPrimaryFieldAssociation() const

Specify the primary field to operate on.

In the cross product operation A x B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetUseCoordinateSystemAsPrimaryField(bool flag)

Specify the primary field to operate on.

In the cross product operation A x B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline bool GetUseCoordinateSystemAsPrimaryField() const

Specify the primary field to operate on.

In the cross product operation A x B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetPrimaryCoordinateSystem(vtkm::Id index)

Specify the primary field to operate on.

In the cross product operation A x B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

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

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline const std::string &GetSecondaryFieldName() const

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::cont::Field::Association GetSecondaryFieldAssociation() const

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetUseCoordinateSystemAsSecondaryField(bool flag)

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline bool GetUseCoordinateSystemAsSecondaryField() const

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetSecondaryCoordinateSystem(vtkm::Id index)

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::Id GetSecondaryCoordinateSystemIndex() const

Specify the secondary field to operate on.

In the cross product operation A x B, B is the secondary field.

The secondary field is an alias for the active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

2.7.13.2. Dot Product

The vtkm::filter::vector_analysis::DotProduct filter computes the dot product of two vector fields for every element in the input data set. The dot product filter computes (PrimaryField . SecondaryField). The dot product computation works for either point or cell centered vector fields.

class DotProduct : public vtkm::filter::Filter

Compute the dot product of vector fields.

The left part of the operand is the “primary” field and the right part of the operand is the “secondary” field (although the dot product is commutative, so the order of primary and secondary seldom matters).

The dot product can operate on vectors of any length.

Public Functions

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

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline const std::string &GetPrimaryFieldName() const

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::cont::Field::Association GetPrimaryFieldAssociation() const

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetUseCoordinateSystemAsPrimaryField(bool flag)

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline bool GetUseCoordinateSystemAsPrimaryField() const

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetPrimaryCoordinateSystem(vtkm::Id coord_idx)

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::Id GetPrimaryCoordinateSystemIndex() const

Specify the primary field to operate on.

In the dot product operation A . B, A is the primary field.

The primary field is an alias for active field index 0. As with any active field, it can be set as a named field or as a coordinate system.

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

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline const std::string &GetSecondaryFieldName() const

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::cont::Field::Association GetSecondaryFieldAssociation() const

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetUseCoordinateSystemAsSecondaryField(bool flag)

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline bool GetUseCoordinateSystemAsSecondaryField() const

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline void SetSecondaryCoordinateSystem(vtkm::Id index)

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

inline vtkm::Id GetSecondaryCoordinateSystemIndex() const

Specify the secondary field to operate on.

In the dot product operation A . B, B is the secondary field.

The secondary field is an alias for active field index 1. As with any active field, it can be set as a named field or as a coordinate system.

2.7.13.3. Gradients

The vtkm::filter::vector_analysis::Gradient filter estimates the gradient of a point based input field for every element in the input data set. The gradient computation can either generate cell center based gradients, which are fast but less accurate, or more accurate but slower point based gradients. The default for the filter is output as cell centered gradients, but can be changed by using the vtkm::filter::vector_analysis::Gradient::SetComputePointGradient() method. The default name for the output fields is “Gradients”, but that can be overridden as always using the vtkm::filter::vector_analysis::Gradient::SetOutputFieldName() method.

class Gradient : public vtkm::filter::Filter

A general filter for gradient estimation.

Estimates the gradient of a point field in a data set. The created gradient array can be determined at either each point location or at the center of each cell.

The default for the filter is output as cell centered gradients. To enable point based gradient computation enable SetComputePointGradient()

If no explicit name for the output field is provided the filter will default to “Gradients”

Public Functions

inline void SetComputePointGradient(bool enable)

Specify whether to compute gradients.

When this flag is on (default is off), the gradient filter will provide a point based gradients, which are significantly more costly since for each point we need to compute the gradient of each cell that uses it.

inline bool GetComputePointGradient() const

Specify whether to compute gradients.

When this flag is on (default is off), the gradient filter will provide a point based gradients, which are significantly more costly since for each point we need to compute the gradient of each cell that uses it.

inline void SetComputeDivergence(bool enable)

Add divergence field to the output data.

The input array must have 3 components to compute this. The default is off.

inline bool GetComputeDivergence() const

Add divergence field to the output data.

The input array must have 3 components to compute this. The default is off.

inline void SetDivergenceName(const std::string &name)

When SetComputeDivergence() is enabled, the result is stored in a field of this name.

If not specified, the name of the field will be Divergence.

inline const std::string &GetDivergenceName() const

When SetComputeDivergence() is enabled, the result is stored in a field of this name.

If not specified, the name of the field will be Divergence.

inline void SetComputeVorticity(bool enable)

Add voriticity/curl field to the output data.

The input array must have 3 components to compute this. The default is off.

inline bool GetComputeVorticity() const

Add voriticity/curl field to the output data.

The input array must have 3 components to compute this. The default is off.

inline void SetVorticityName(const std::string &name)

When SetComputeVorticity() is enabled, the result is stored in a field of this name.

If not specified, the name of the field will be Vorticity.

inline const std::string &GetVorticityName() const

When SetComputeVorticity() is enabled, the result is stored in a field of this name.

If not specified, the name of the field will be Vorticity.

inline void SetComputeQCriterion(bool enable)

Add Q-criterion field to the output data.

The input array must have 3 components to compute this. The default is off.

inline bool GetComputeQCriterion() const

Add Q-criterion field to the output data.

The input array must have 3 components to compute this. The default is off.

inline void SetQCriterionName(const std::string &name)

When SetComputeQCriterion() is enabled, the result is stored in a field of this name.

If not specified, the name of the field will be QCriterion.

inline const std::string &GetQCriterionName() const

When SetComputeQCriterion() is enabled, the result is stored in a field of this name.

If not specified, the name of the field will be QCriterion.

inline void SetComputeGradient(bool enable)

Add gradient field to the output data.

The name of the array will be Gradients unless otherwise specified with SetOutputFieldName and will be a cell field unless ComputePointGradient() is enabled. It is useful to turn this off when you are only interested in the results of Divergence, Vorticity, or QCriterion. The default is on.

inline bool GetComputeGradient() const

Add gradient field to the output data.

The name of the array will be Gradients unless otherwise specified with SetOutputFieldName and will be a cell field unless ComputePointGradient() is enabled. It is useful to turn this off when you are only interested in the results of Divergence, Vorticity, or QCriterion. The default is on.

inline void SetColumnMajorOrdering()

Make the vector gradient output format be in FORTRAN Column-major order.

This is only used when the input field is a vector field. Enabling column-major is important if integrating with other projects such as VTK. Default: Row Order.

inline void SetRowMajorOrdering()

Make the vector gradient output format be in C Row-major order.

This is only used when the input field is a vector field. Default: Row Order.

2.7.13.4. Surface Normals

The vtkm::filter::vector_analysis::SurfaceNormals filter computes the surface normals of a polygonal data set at its points and/or cells. The filter takes a data set as input and by default, uses the active coordinate system to compute the normals.

class SurfaceNormals : public vtkm::filter::Filter

Computes normals for polygonal mesh.

This filter computes surface normals on points and/or cells of a polygonal dataset. The cell normals are faceted and are computed based on the plane where a face lies. The point normals are smooth normals, computed by averaging the face normals of incident cells. The normals will be consistently oriented to point in the direction of the same connected surface if possible.

The point and cell normals may be oriented to a point outside of the manifold surface by turning on the auto orient normals option (SetAutoOrientNormals()), or they may point inward by also setting flip normals (SetFlipNormals()) to true.

Triangle vertices will be reordered to be wound counter-clockwise around the cell normals when the consistency option (SetConsistency()) is enabled.

For non-polygonal cells, a zeroed vector is assigned. The point normals are computed by averaging the cell normals of the incident cells of each point.

The default name for the output fields is Normals, but that can be overridden using the SetCellNormalsName() and SetPointNormalsName() methods. The filter will also respect the name in SetOutputFieldName() if neither of the others are set.

Public Functions

SurfaceNormals()

Create SurfaceNormals filter.

This calls this->SetUseCoordinateSystemAsField(true) since that is the most common use-case for surface normals.

inline void SetGenerateCellNormals(bool value)

Specify whether cell normals should be generated.

Default is off.

inline bool GetGenerateCellNormals() const

Specify whether cell normals should be generated.

Default is off.

inline void SetNormalizeCellNormals(bool value)

Specify whether the cell normals should be normalized.

Default value is true. The intended use case of this flag is for faster, approximate point normals generation by skipping the normalization of the face normals. Note that when set to false, the result cell normals will not be unit length normals and the point normals will be different.

inline bool GetNormalizeCellNormals() const

Specify whether the cell normals should be normalized.

Default value is true. The intended use case of this flag is for faster, approximate point normals generation by skipping the normalization of the face normals. Note that when set to false, the result cell normals will not be unit length normals and the point normals will be different.

inline void SetGeneratePointNormals(bool value)

Specify whether the point normals should be generated.

Default is on.

inline bool GetGeneratePointNormals() const

Specify whether the point normals should be generated.

Default is on.

inline void SetCellNormalsName(const std::string &name)

Specify the name of the cell normals field.

Default is Normals.

inline const std::string &GetCellNormalsName() const

Specify the name of the cell normals field.

Default is Normals.

inline void SetPointNormalsName(const std::string &name)

Specify the name of the point normals field.

Default is Normals.

inline const std::string &GetPointNormalsName() const

Specify the name of the point normals field.

Default is Normals.

inline void SetAutoOrientNormals(bool v)

Specify whether to orient the normals outwards from the surface.

This requires a closed manifold surface or the behavior is undefined. This option is expensive but might be necessary for rendering. To make the normals point inward, set FlipNormals to true. Default is off.

inline bool GetAutoOrientNormals() const

Specify whether to orient the normals outwards from the surface.

This requires a closed manifold surface or the behavior is undefined. This option is expensive but might be necessary for rendering. To make the normals point inward, set FlipNormals to true. Default is off.

inline void SetFlipNormals(bool v)

Specify the direction to point normals when SetAutoOrientNormals() is true.

When this flag is false (the default), the normals will be oriented to point outward. When the flag is true, the normals will point inward. This option has no effect if auto orient normals is off.

inline bool GetFlipNormals() const

Specify the direction to point normals when SetAutoOrientNormals() is true.

When this flag is false (the default), the normals will be oriented to point outward. When the flag is true, the normals will point inward. This option has no effect if auto orient normals is off.

inline void SetConsistency(bool v)

Specify whtehr polygon winding should be made consistent with normal orientation.

Triangles are wound such that their points are counter-clockwise around the generated cell normal. Default is true. This currently only affects triangles. This is only applied when cell normals are generated.

inline bool GetConsistency() const

Specify whtehr polygon winding should be made consistent with normal orientation.

Triangles are wound such that their points are counter-clockwise around the generated cell normal. Default is true. This currently only affects triangles. This is only applied when cell normals are generated.

2.7.13.5. Vector Magnitude

The vtkm::filter::vector_analysis::VectorMagnitude filter takes a field comprising vectors and computes the magnitude for each vector. The vector field is selected as usual with the vtkm::filter::vector_analysis::VectorMagnitude::SetActiveField() method. The default name for the output field is magnitude, but that can be overridden as always using the vtkm::filter::vector_analysis::VectorMagnitude::SetOutputFieldName() method.

class VectorMagnitude : public vtkm::filter::Filter

Compute the magnitudes of a vector field.

The vector field is selected with the SetActiveField() method. The default name for the output field is magnitude, but that can be overridden using the SetOutputFieldName() method.

2.7.14. ZFP Compression

vtkm::filter::zfp::ZFPCompressor1D, vtkm::filter::zfp::ZFPCompressor2D, and vtkm::filter::zfp::ZFPCompressor3D are a set of filters that take a 1D, 2D, and 3D field, respectively, and compresses the values using the compression algorithm ZFP.

The field is selected as usual with the vtkm::filter::zfp::ZFPCompressor3D::SetActiveField() method. The rate of compression is set using vtkm::filter::zfp::ZFPCompressor3D::SetRate(). The default name for the output field is compressed.

class ZFPCompressor1D : public vtkm::filter::Filter

Compress a scalar field using ZFP.

Takes as input a 1D array and generates an output of compressed data.

Warning

This filter currently only supports 1D structured cell sets.

Public Functions

inline void SetRate(vtkm::Float64 _rate)

Specifies the rate of compression.

inline vtkm::Float64 GetRate()

Specifies the rate of compression.

class ZFPCompressor2D : public vtkm::filter::Filter

Compress a scalar field using ZFP.

Takes as input a 2D array and generates an output of compressed data.

Warning

This filter is currently only supports 2D structured cell sets.

Public Functions

inline void SetRate(vtkm::Float64 _rate)

Specifies the rate of compression.

inline vtkm::Float64 GetRate()

Specifies the rate of compression.

class ZFPCompressor3D : public vtkm::filter::Filter

Compress a scalar field using ZFP.

Takes as input a 3D array and generates an output of compressed data.

Warning

This filter is currently only supports 3D structured cell sets.

Public Functions

inline void SetRate(vtkm::Float64 _rate)

Specifies the rate of compression.

inline vtkm::Float64 GetRate()

Specifies the rate of compression.

vtkm::filter::zfp::ZFPDecompressor1D, vtkm::filter::zfp::ZFPDecompressor2D, and vtkm::filter::zfp::ZFPDecompressor3D are a set of filters that take a compressed 1D, 2D, and 3D field, respectively, and decompress the values using the compression algorithm ZFP.

The field is selected as usual with the vtkm::filter::zfp::ZFPDecompressor3D::SetActiveField() method. The rate of compression is set using vtkm::filter::zfp::ZFPDecompressor3D::SetRate(). The default name for the output field is decompressed.

class ZFPDecompressor1D : public vtkm::filter::Filter

Decompress a scalar field using ZFP.

Takes as input a 1D compressed array and generates the decompressed version of the data.

Warning

This filter is currently only supports 1D structured cell sets.

Public Functions

inline void SetRate(vtkm::Float64 _rate)

Specifies the rate of compression.

inline vtkm::Float64 GetRate()

Specifies the rate of compression.

class ZFPDecompressor2D : public vtkm::filter::Filter

Decompress a scalar field using ZFP.

Takes as input a 2D compressed array and generates the decompressed version of the data.

Warning

This filter is currently only supports 2D structured cell sets.

Public Functions

inline void SetRate(vtkm::Float64 _rate)

Specifies the rate of compression.

inline vtkm::Float64 GetRate()

Specifies the rate of compression.

class ZFPDecompressor3D : public vtkm::filter::Filter

Decompress a scalar field using ZFP.

Takes as input a 3D compressed array and generates the decompressed version of the data.

Warning

This filter is currently only supports 3D structured cell sets.

Public Functions

inline void SetRate(vtkm::Float64 _rate)

Specifies the rate of compression.

inline vtkm::Float64 GetRate()

Specifies the rate of compression.