2.4. Data Sets

A data set, implemented with the vtkm::cont::DataSet class, contains and manages the geometric data structures that VTK‑m operates on.

class DataSet

Contains and manages the geometric data structures that VTK-m operates on.

A DataSet is the main data structure used by VTK-m to pass data in and out of filters, rendering, and other components. A data set comprises the following 3 data structures.

  • CellSet A cell set describes topological connections. A cell set defines some number of points in space and how they connect to form cells, filled regions of space. A data set has exactly one cell set.

  • Field A field describes numerical data associated with the topological elements in a cell set. The field is represented as an array, and each entry in the field array corresponds to a topological element (point, edge, face, or cell). Together the cell set topology and discrete data values in the field provide an interpolated function throughout the volume of space covered by the data set. A cell set can have any number of fields.

  • CoordinateSystem A coordinate system is a special field that describes the physical location of the points in a data set. Although it is most common for a data set to contain a single coordinate system, VTK-m supports data sets with no coordinate system such as abstract data structures like graphs that might not have positions in a space. DataSet also supports multiple coordinate systems for data that have multiple representations for position. For example, geospatial data could simultaneously have coordinate systems defined by 3D position, latitude-longitude, and any number of 2D projections.

In addition to the base vtkm::cont::DataSet, VTK‑m provides vtkm::cont::PartitionedDataSet to represent data partitioned into multiple domains. A vtkm::cont::PartitionedDataSet is implemented as a collection of vtkm::cont::DataSet objects. Partitioned data sets are described later in Section 2.4.5 (Partitioned Data Sets).

2.4.1. Building Data Sets

Before we go into detail on the cell sets, fields, and coordinate systems that make up a data set in VTK‑m, let us first discuss how to build a data set. One simple way to build a data set is to load data from a file using the vtkm::io module. Reading files is discussed in detail in Chapter 2.5 (File I/O).

This section describes building data sets of different types using a set of classes named DataSetBuilder*, which provide a convenience layer on top of vtkm::cont::DataSet to make it easier to create data sets.

Did You Know?

To simplify the introduction of vtkm::cont::DataSet objects, this section uses the simplest mechanisms. In many cases this involves loading data in a std::vector and passing that to VTK‑m, which usually causes the data to be copied. This is not the most efficient method to load data into VTK‑m. Although it is sufficient for small data or data that come from a “slow” source, such as a file, it might be a bottleneck for large data generated by another library. It is possible to adapt VTK‑m’s vtkm::cont::DataSet to externally defined data. This is done by wrapping existing data into what is called ArrayHandle, but this is a more advanced topic that will not be addressed in this chapter. ArrayHandle objects are introduced in Chapter 3.2 (Basic Array Handles) and more adaptive techniques are described in later chapters.

2.4.1.1. Creating Uniform Grids

Uniform grids are meshes that have a regular array structure with points uniformly spaced parallel to the axes. Uniform grids are also sometimes called regular grids or images.

The vtkm::cont::DataSetBuilderUniform class can be used to easily create 2- or 3-dimensional uniform grids. vtkm::cont::DataSetBuilderUniform has several versions of a method named vtkm::cont::DataSetBuilderUniform::Create() that takes the number of points in each dimension, the origin, and the spacing. The origin is the location of the first point of the data (in the lower left corner), and the spacing is the distance between points in the x, y, and z directions.

class DataSetBuilderUniform

Public Static Functions

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::Id &dimension, const T &origin, const T &spacing, const std::string &coordNm = "coords")

Create a 1D uniform DataSet.

Parameters:
  • dimension[in] The size of the grid. The dimensions are specified based on the number of points (as opposed to the number of cells).

  • origin[in] The origin of the data. This is the point coordinate with the minimum value in all dimensions.

  • spacing[in] The uniform distance between adjacent points.

  • coordNm[in] (optional) The name to register the coordinates as.

static vtkm::cont::DataSet Create(const vtkm::Id &dimension, const std::string &coordNm = "coords")

Create a 1D uniform DataSet.

The origin is set to 0 and the spacing is set to 1.

Parameters:
  • dimension[in] The size of the grid. The dimensions are specified based on the number of points (as opposed to the number of cells).

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::Id2 &dimensions, const vtkm::Vec<T, 2> &origin, const vtkm::Vec<T, 2> &spacing, const std::string &coordNm = "coords")

Create a 2D uniform DataSet.

Parameters:
  • dimensions[in] The size of the grid. The dimensions are specified based on the number of points (as opposed to the number of cells).

  • origin[in] The origin of the data. This is the point coordinate with the minimum value in all dimensions.

  • spacing[in] The uniform distance between adjacent points.

  • coordNm[in] (optional) The name to register the coordinates as.

static vtkm::cont::DataSet Create(const vtkm::Id2 &dimensions, const std::string &coordNm = "coords")

Create a 2D uniform DataSet.

The origin is set to (0,0) and the spacing is set to (1,1).

Parameters:
  • dimensions[in] The size of the grid. The dimensions are specified based on the number of points (as opposed to the number of cells).

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::Id3 &dimensions, const vtkm::Vec<T, 3> &origin, const vtkm::Vec<T, 3> &spacing, const std::string &coordNm = "coords")

Create a 3D uniform DataSet.

Parameters:
  • dimensions[in] The size of the grid. The dimensions are specified based on the number of points (as opposed to the number of cells).

  • origin[in] The origin of the data. This is the point coordinate with the minimum value in all dimensions.

  • spacing[in] The uniform distance between adjacent points.

  • coordNm[in] (optional) The name to register the coordinates as.

static vtkm::cont::DataSet Create(const vtkm::Id3 &dimensions, const std::string &coordNm = "coords")

Create a 3D uniform DataSet.

The origin is set to (0,0,0) and the spacing is set to (1,1,1).

Parameters:
  • dimensions[in] The size of the grid. The dimensions are specified based on the number of points (as opposed to the number of cells).

  • coordNm[in] (optional) The name to register the coordinates as.

The following example creates a vtkm::cont::DataSet containing a uniform grid of \(101 \times 101 \times 26\) points.

Example 2.3 Creating a uniform grid.}{.cxx}
1  vtkm::cont::DataSetBuilderUniform dataSetBuilder;
2
3  vtkm::cont::DataSet dataSet = dataSetBuilder.Create(vtkm::Id3(101, 101, 26));

If not specified, the origin will be at the coordinates \((0,0,0)\) and the spacing will be \(1\) in each direction. Thus, in the previous example the width, height, and depth of the mesh in physical space will be \(100\), \(100\), and :math`25`, respectively, and the mesh will be centered at \((50, 50, 12.5)\). Let us say we actually want a mesh of the same dimensions, but we want the \(z\) direction to be stretched out so that the mesh will be the same size in each direction, and we want the mesh centered at the origin.

Example 2.4 Creating a uniform grid with custom origin and spacing.
1  vtkm::cont::DataSetBuilderUniform dataSetBuilder;
2
3  vtkm::cont::DataSet dataSet = dataSetBuilder.Create(vtkm::Id3(101, 101, 26),
4                                                      vtkm::Vec3f(-50.0, -50.0, -50.0),
5                                                      vtkm::Vec3f(1.0, 1.0, 4.0));

2.4.1.2. Creating Rectilinear Grids

A rectilinear grid is similar to a uniform grid except that a rectilinear grid can adjust the spacing between adjacent grid points. This allows the rectilinear grid to have tighter sampling in some areas of space, but the points are still constrained to be aligned with the axes and each other. The irregular spacing of a rectilinear grid is specified by providing a separate array each for the x, y, and z coordinates.

The vtkm::cont::DataSetBuilderRectilinear class can be used to easily create 2- or 3-dimensional rectilinear grids. vtkm::cont::DataSetBuilderRectilinear has several versions of a method named vtkm::cont::DataSetBuilderRectilinear::Create() that takes these coordinate arrays and builds a vtkm::cont::DataSet out of them. The arrays can be supplied as either standard C arrays or as std::vector objects, in which case the data in the arrays are copied into the vtkm::cont::DataSet. These arrays can also be passed as vtkm::cont::ArrayHandle objects (introduced later in this book), in which case the data are shallow copied.

class DataSetBuilderRectilinear

Public Static Functions

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<T> &xvals, const std::string &coordNm = "coords")

Create a 1D retilinear DataSet.

A rectilinear grid is specified with a scalar array for the point coordinates in the x direction. In this form, the coordinate array is specified with std::vector. The data is copied from the std::vector.

Parameters:
  • xvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(vtkm::Id nx, T *xvals, const std::string &coordNm = "coords")

Create a 1D retilinear DataSet.

A rectilinear grid is specified with a scalar array for the point coordinates in the x direction. In this form, the coordinate array is specified with a standard C array. The data is copied from the array.

Parameters:
  • nx[in] The size of the grid in the x direction (and length of the xvals array).

  • xvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle<T> &xvals, const std::string &coordNm = "coords")

Create a 1D retilinear DataSet.

A rectilinear grid is specified with a scalar array for the point coordinates in the x direction. In this form, the coordinate array is specified with vtkm::cont::ArrayHandle. The ArrayHandle is shared with the DataSet, so changing the ArrayHandle changes the DataSet.

Parameters:
  • xvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<T> &xvals, const std::vector<T> &yvals, const std::string &coordNm = "coords")

Create a 2D retilinear DataSet.

A rectilinear grid is specified with separate arrays for the point coordinates in the x and y directions. In this form, the coordinate arrays are specified with std::vector. The data is copied from the std::vectors.

Parameters:
  • xvals[in] An array of coordinates to use along the x dimension.

  • yvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(vtkm::Id nx, vtkm::Id ny, T *xvals, T *yvals, const std::string &coordNm = "coords")

Create a 2D retilinear DataSet.

A rectilinear grid is specified with separate arrays for the point coordinates in the x and y directions. In this form, the coordinate arrays are specified with standard C arrays. The data is copied from the arrays.

Parameters:
  • nx[in] The size of the grid in the x direction (and length of the xvals array).

  • ny[in] The size of the grid in the x direction (and length of the yvals array).

  • xvals[in] An array of coordinates to use along the x dimension.

  • yvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle<T> &xvals, const vtkm::cont::ArrayHandle<T> &yvals, const std::string &coordNm = "coords")

Create a 2D retilinear DataSet.

A rectilinear grid is specified with separate arrays for the point coordinates in the x and y directions. In this form, the coordinate arrays are specified with vtkm::cont::ArrayHandle. The ArrayHandles are shared with the DataSet, so changing the ArrayHandles changes the DataSet.

Parameters:
  • xvals[in] An array of coordinates to use along the x dimension.

  • yvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(vtkm::Id nx, vtkm::Id ny, vtkm::Id nz, T *xvals, T *yvals, T *zvals, const std::string &coordNm = "coords")

Create a 3D retilinear DataSet.

A rectilinear grid is specified with separate arrays for the point coordinates in the x, y, and z directions. In this form, the coordinate arrays are specified with standard C arrays. The data is copied from the arrays.

Parameters:
  • nx[in] The size of the grid in the x direction (and length of the xvals array).

  • ny[in] The size of the grid in the x direction (and length of the yvals array).

  • nz[in] The size of the grid in the x direction (and length of the zvals array).

  • xvals[in] An array of coordinates to use along the x dimension.

  • yvals[in] An array of coordinates to use along the x dimension.

  • zvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<T> &xvals, const std::vector<T> &yvals, const std::vector<T> &zvals, const std::string &coordNm = "coords")

Create a 3D retilinear DataSet.

A rectilinear grid is specified with separate arrays for the point coordinates in the x, y, and z directions. In this form, the coordinate arrays are specified with std::vector. The data is copied from the std::vectors.

Parameters:
  • xvals[in] An array of coordinates to use along the x dimension.

  • yvals[in] An array of coordinates to use along the x dimension.

  • zvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle<T> &xvals, const vtkm::cont::ArrayHandle<T> &yvals, const vtkm::cont::ArrayHandle<T> &zvals, const std::string &coordNm = "coords")

Create a 3D retilinear DataSet.

A rectilinear grid is specified with separate arrays for the point coordinates in the x, y, and z directions. In this form, the coordinate arrays are specified with vtkm::cont::ArrayHandle. The ArrayHandles are shared with the DataSet, so changing the ArrayHandles changes the DataSet.

Parameters:
  • xvals[in] An array of coordinates to use along the x dimension.

  • yvals[in] An array of coordinates to use along the x dimension.

  • zvals[in] An array of coordinates to use along the x dimension.

  • coordNm[in] (optional) The name to register the coordinates as.

The following example creates a vtkm::cont::DataSet containing a rectilinear grid with \(201 \times 201 \times 101\) points with different irregular spacing along each axis.

Example 2.5 Creating a rectilinear grid.
 1  // Make x coordinates range from -4 to 4 with tighter spacing near 0.
 2  std::vector<vtkm::Float32> xCoordinates;
 3  for (vtkm::Float32 x = -2.0f; x <= 2.0f; x += 0.02f)
 4  {
 5    xCoordinates.push_back(vtkm::CopySign(x * x, x));
 6  }
 7
 8  // Make y coordinates range from 0 to 2 with tighter spacing near 2.
 9  std::vector<vtkm::Float32> yCoordinates;
10  for (vtkm::Float32 y = 0.0f; y <= 4.0f; y += 0.02f)
11  {
12    yCoordinates.push_back(vtkm::Sqrt(y));
13  }
14
15  // Make z coordinates rangefrom -1 to 1 with even spacing.
16  std::vector<vtkm::Float32> zCoordinates;
17  for (vtkm::Float32 z = -1.0f; z <= 1.0f; z += 0.02f)
18  {
19    zCoordinates.push_back(z);
20  }
21
22  vtkm::cont::DataSetBuilderRectilinear dataSetBuilder;
23
24  vtkm::cont::DataSet dataSet =
25    dataSetBuilder.Create(xCoordinates, yCoordinates, zCoordinates);

2.4.1.3. Creating Explicit Meshes

An explicit mesh is an arbitrary collection of cells with arbitrary connections. It can have multiple different types of cells. Explicit meshes are also known as unstructured grids. Explicit meshes can contain cells of different shapes. The shapes that VTK‑m currently supports are listed in Figure 2.1. Each shape is identified using either a numeric identifier, provided by VTK‑m with identifiers of the form vtkm::CELL_SHAPE_* or special tag structures of the form vtkm::CellSetTag*. Cell shapes are discussed in detail in Chapter 4.7 (Working with Cells).

_images/CellConnections.png

Figure 2.1 Basic Cell Shapes.

_images/ExplicitCellConnections.png

Figure 2.2 An example explicit mesh.

The cells of an explicit mesh are defined with the following 3 arrays, which are depicted graphically in Figure 2.2.

Shapes

An array of ids identifying the shape of the cell. Each value is a vtkm::UInt8 and should be set to one of the vtkm::CELL_SHAPE_* constants. The shapes and their identifiers are shown in Figure 2.1. The size of this array is equal to the number of cells in the set.

Connectivity

An array that lists all the points that comprise each cell. Each entry in the array is a vtkm::Id giving the point id associated with a vertex of a cell. The points for each cell are given in a prescribed order for each shape, which is also shown in Figure 2.1. The point indices are stored consecutively from the first cell to the last.

Offsets

An array of vtkm::Id’s pointing to the index in the connectivity array where the points for a particular cell starts. The size of this array is equal to the number of cells in the set plus 1. The first entry is expected to be 0 (since the connectivity of the first cell is at the start of the connectivity array). The last entry, which does not correspond to any cell, should be the size of the connectivity array.

One important item that is missing from this list of arrays is a count of the number of indices associated with each cell. This is not explicitly represented in VTK‑m’s mesh structure because it can be implicitly derived from the offsets array by subtracting consecutive entries. However, it is usually the case when building an explicit mesh that you will have an array of these counts rather than the offsets. It is for this reason that VTK‑m contains mechanisms to build an explicit data set with a “num indices” arrays rather than an offsets array.

The vtkm::cont::DataSetBuilderExplicit class can be used to create data sets with explicit meshes. vtkm::cont::DataSetBuilderExplicit has several versions of a method named vtkm::cont::DataSetBuilderExplicit::Create(). Generally, these methods take the shapes, number of indices, and connectivity arrays as well as an array of point coordinates.

class DataSetBuilderExplicit

Public Static Functions

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<T> &xVals, const std::vector<vtkm::UInt8> &shapes, const std::vector<vtkm::IdComponent> &numIndices, const std::vector<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 1D DataSet with arbitrary cell connectivity.

The cell connectivity is specified with arrays defining the shape and point connections of each cell. In this form, the cell connectivity and coordinates are specified as std::vector and the data will be copied to create the data object.

Parameters:
  • xVals[in] An array providing the x coordinate of each point.

  • shapes[in] An array of shapes for each cell. Each entry should be one of the vtkm::CELL_SHAPE_* values identifying the shape of the corresponding cell.

  • numIndices[in] An array containing for each cell the number of points incident on that cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by the numIndices array. These variable length arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<T> &xVals, const std::vector<T> &yVals, const std::vector<vtkm::UInt8> &shapes, const std::vector<vtkm::IdComponent> &numIndices, const std::vector<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 2D DataSet with arbitrary cell connectivity.

The cell connectivity is specified with arrays defining the shape and point connections of each cell. In this form, the cell connectivity and coordinates are specified as std::vector and the data will be copied to create the data object.

Parameters:
  • xVals[in] An array providing the x coordinate of each point.

  • yVals[in] An array providing the x coordinate of each point.

  • shapes[in] An array of shapes for each cell. Each entry should be one of the vtkm::CELL_SHAPE_* values identifying the shape of the corresponding cell.

  • numIndices[in] An array containing for each cell the number of points incident on that cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by the numIndices array. These variable length arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<T> &xVals, const std::vector<T> &yVals, const std::vector<T> &zVals, const std::vector<vtkm::UInt8> &shapes, const std::vector<vtkm::IdComponent> &numIndices, const std::vector<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 3D DataSet with arbitrary cell connectivity.

The cell connectivity is specified with arrays defining the shape and point connections of each cell. In this form, the cell connectivity and coordinates are specified as std::vector and the data will be copied to create the data object.

Parameters:
  • xVals[in] An array providing the x coordinate of each point.

  • yVals[in] An array providing the x coordinate of each point.

  • zVals[in] An array providing the x coordinate of each point.

  • shapes[in] An array of shapes for each cell. Each entry should be one of the vtkm::CELL_SHAPE_* values identifying the shape of the corresponding cell.

  • numIndices[in] An array containing for each cell the number of points incident on that cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by the numIndices array. These variable length arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const std::vector<vtkm::Vec<T, 3>> &coords, const std::vector<vtkm::UInt8> &shapes, const std::vector<vtkm::IdComponent> &numIndices, const std::vector<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 3D DataSet with arbitrary cell connectivity.

The cell connectivity is specified with arrays defining the shape and point connections of each cell. In this form, the cell connectivity and coordinates are specified as std::vector and the data will be copied to create the data object.

Parameters:
  • coords[in] An array of point coordinates.

  • shapes[in] An array of shapes for each cell. Each entry should be one of the vtkm::CELL_SHAPE_* values identifying the shape of the corresponding cell.

  • numIndices[in] An array containing for each cell the number of points incident on that cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by the numIndices array. These variable length arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

template<typename T>
static inline vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle<vtkm::Vec<T, 3>> &coords, const vtkm::cont::ArrayHandle<vtkm::UInt8> &shapes, const vtkm::cont::ArrayHandle<vtkm::IdComponent> &numIndices, const vtkm::cont::ArrayHandle<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 3D DataSet with arbitrary cell connectivity.

The cell connectivity is specified with arrays defining the shape and point connections of each cell. In this form, the cell connectivity and coordinates are specified as ArrayHandle and the memory will be shared with the created data object. That said, the DataSet construction will generate a new array for offsets.

Parameters:
  • coords[in] An array of point coordinates.

  • shapes[in] An array of shapes for each cell. Each entry should be one of the vtkm::CELL_SHAPE_* values identifying the shape of the corresponding cell.

  • numIndices[in] An array containing for each cell the number of points incident on that cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by the numIndices array. These variable length arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

template<typename T, typename CellShapeTag>
static inline vtkm::cont::DataSet Create(const std::vector<vtkm::Vec<T, 3>> &coords, CellShapeTag tag, vtkm::IdComponent numberOfPointsPerCell, const std::vector<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 3D DataSet with arbitrary cell connectivity for a single cell type.

The cell connectivity is specified with an array defining the point connections of each cell. All the cells in the DataSet are of the same shape and contain the same number of incident points. In this form, the cell connectivity and coordinates are specified as std::vector and the data will be copied to create the data object.

Parameters:
  • coords[in] An array of point coordinates.

  • tag[in] A tag object representing the shape of all the cells in the mesh. Cell shape tag objects have a name of the form vtkm::CellShapeTag* such as vtkm::CellShapeTagTriangle or vtkm::CellShapeTagHexahedron. To specify a cell shape determined at runtime, use vtkm::CellShapeTagGeneric.

  • numberOfPointsPerCell[in] The number of points that are incident to each cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by numberOfPointsPerCell. These short arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

template<typename T, typename CellShapeTag>
static inline vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle<vtkm::Vec<T, 3>> &coords, CellShapeTag tag, vtkm::IdComponent numberOfPointsPerCell, const vtkm::cont::ArrayHandle<vtkm::Id> &connectivity, const std::string &coordsNm = "coords")

Create a 3D DataSet with arbitrary cell connectivity for a single cell type.

The cell connectivity is specified with an array defining the point connections of each cell. All the cells in the DataSet are of the same shape and contain the same number of incident points. In this form, the cell connectivity and coordinates are specified as ArrayHandle and the memory will be shared with the created data object.

Parameters:
  • coords[in] An array of point coordinates.

  • tag[in] A tag object representing the shape of all the cells in the mesh. Cell shape tag objects have a name of the form vtkm::CellShapeTag* such as vtkm::CellShapeTagTriangle or vtkm::CellShapeTagHexahedron. To specify a cell shape determined at runtime, use vtkm::CellShapeTagGeneric.

  • numberOfPointsPerCell[in] The number of points that are incident to each cell.

  • connectivity[in] An array specifying for each cell the indicies of points incident on each cell. Each cell has a short array of indices that reference points in the coords array. The length of each of these short arrays is specified by numberOfPointsPerCell. These short arrays are tightly packed together in this connectivity array.

  • coordsNm[in] (optional) The name to register the coordinates as.

The following example creates a mesh like the one shown in Figure 2.2.

Example 2.6 Creating an explicit mesh with vtkm::cont::DataSetBuilderExplicit.
 1  // Array of point coordinates.
 2  std::vector<vtkm::Vec3f_32> pointCoordinates;
 3  pointCoordinates.push_back(vtkm::Vec3f_32(1.1f, 0.0f, 0.0f));
 4  pointCoordinates.push_back(vtkm::Vec3f_32(0.2f, 0.4f, 0.0f));
 5  pointCoordinates.push_back(vtkm::Vec3f_32(0.9f, 0.6f, 0.0f));
 6  pointCoordinates.push_back(vtkm::Vec3f_32(1.4f, 0.5f, 0.0f));
 7  pointCoordinates.push_back(vtkm::Vec3f_32(1.8f, 0.3f, 0.0f));
 8  pointCoordinates.push_back(vtkm::Vec3f_32(0.4f, 1.0f, 0.0f));
 9  pointCoordinates.push_back(vtkm::Vec3f_32(1.0f, 1.2f, 0.0f));
10  pointCoordinates.push_back(vtkm::Vec3f_32(1.5f, 0.9f, 0.0f));
11
12  // Array of shapes.
13  std::vector<vtkm::UInt8> shapes;
14  shapes.push_back(vtkm::CELL_SHAPE_TRIANGLE);
15  shapes.push_back(vtkm::CELL_SHAPE_QUAD);
16  shapes.push_back(vtkm::CELL_SHAPE_TRIANGLE);
17  shapes.push_back(vtkm::CELL_SHAPE_POLYGON);
18  shapes.push_back(vtkm::CELL_SHAPE_TRIANGLE);
19
20  // Array of number of indices per cell.
21  std::vector<vtkm::IdComponent> numIndices;
22  numIndices.push_back(3);
23  numIndices.push_back(4);
24  numIndices.push_back(3);
25  numIndices.push_back(5);
26  numIndices.push_back(3);
27
28  // Connectivity array.
29  std::vector<vtkm::Id> connectivity;
30  connectivity.push_back(0); // Cell 0
31  connectivity.push_back(2);
32  connectivity.push_back(1);
33  connectivity.push_back(0); // Cell 1
34  connectivity.push_back(4);
35  connectivity.push_back(3);
36  connectivity.push_back(2);
37  connectivity.push_back(1); // Cell 2
38  connectivity.push_back(2);
39  connectivity.push_back(5);
40  connectivity.push_back(2); // Cell 3
41  connectivity.push_back(3);
42  connectivity.push_back(7);
43  connectivity.push_back(6);
44  connectivity.push_back(5);
45  connectivity.push_back(3); // Cell 4
46  connectivity.push_back(4);
47  connectivity.push_back(7);
48
49  // Copy these arrays into a DataSet.
50  vtkm::cont::DataSetBuilderExplicit dataSetBuilder;
51
52  vtkm::cont::DataSet dataSet =
53    dataSetBuilder.Create(pointCoordinates, shapes, numIndices, connectivity);

Often it is awkward to build your own arrays and then pass them to vtkm::cont::DataSetBuilderExplicit. There also exists an alternate builder class named vtkm::cont::DataSetBuilderExplicitIterative that allows you to specify each cell and point one at a time rather than all at once. This is done by calling one of the versions of vtkm::cont::DataSetBuilderExplicitIterative::AddPoint() and one of the versions of vtkm::cont::DataSetBuilderExplicitIterative::AddCell() for each point and cell, respectively.

class DataSetBuilderExplicitIterative

Helper class to build a DataSet by iteratively adding points and cells.

This class allows you to specify a DataSet by adding points and cells one at a time.

Public Functions

void Begin(const std::string &coordName = "coords")

Begin defining points and cells of a DataSet.

The state of this object is initialized to be ready to use AddPoint and AddCell methods.

Parameters:

coordName[in] (optional) The name to register the coordinates as.

vtkm::Id AddPoint(const vtkm::Vec3f &pt)

Add a point to the DataSet.

Parameters:

pt[in] The coordinates of the point to add.

Returns:

The index of the newly created point.

template<typename T>
inline vtkm::Id AddPoint(const vtkm::Vec<T, 3> &pt)

Add a point to the DataSet.

Parameters:

pt[in] The coordinates of the point to add.

Returns:

The index of the newly created point.

vtkm::Id AddPoint(const vtkm::FloatDefault &x, const vtkm::FloatDefault &y, const vtkm::FloatDefault &z = 0)

Add a point to the DataSet.

Parameters:
  • x[in] The x coordinate of the newly created point.

  • y[in] The y coordinate of the newly created point.

  • z[in] The z coordinate of the newly created point.

Returns:

The index of the newly created point.

template<typename T>
inline vtkm::Id AddPoint(const T &x, const T &y, const T &z = 0)

Add a point to the DataSet.

Parameters:
  • x[in] The x coordinate of the newly created point.

  • y[in] The y coordinate of the newly created point.

  • z[in] The z coordinate of the newly created point.

Returns:

The index of the newly created point.

void AddCell(const vtkm::UInt8 &shape, const std::vector<vtkm::Id> &conn)

Add a cell to the DataSet.

Parameters:
  • shape[in] Identifies the shape of the cell. Use one of the vtkm::CELL_SHAPE_* values.

  • conn[in] List of indices to the incident points.

void AddCell(const vtkm::UInt8 &shape, const vtkm::Id *conn, const vtkm::IdComponent &n)

Add a cell to the DataSet.

Parameters:
  • shape[in] Identifies the shape of the cell. Use one of the vtkm::CELL_SHAPE_* values.

  • conn[in] List of indices to the incident points.

  • n[in] The number of incident points (and the length of the conn array).

void AddCell(vtkm::UInt8 shape)

Start adding a cell to the DataSet.

The incident points are later added one at a time using AddCellPoint. The cell is completed the next time AddCell or Create is called.

Parameters:

shape[in] Identifies the shape of the cell. Use one of the

void AddCellPoint(vtkm::Id pointIndex)

Add an incident point to the current cell.

Parameters:

pointIndex[in] Index to the incident point.

vtkm::cont::DataSet Create()

Produce the DataSet.

The points and cells previously added are finalized and the resulting DataSet is returned.

The next example also builds the mesh shown in Figure 2.2 except this time using vtkm::cont::DataSetBuilderExplicitIterative.

Example 2.7 Creating an explicit mesh with vtkm::cont::DataSetBuilderExplicitIterative.
 1  vtkm::cont::DataSetBuilderExplicitIterative dataSetBuilder;
 2
 3  dataSetBuilder.AddPoint(1.1, 0.0, 0.0);
 4  dataSetBuilder.AddPoint(0.2, 0.4, 0.0);
 5  dataSetBuilder.AddPoint(0.9, 0.6, 0.0);
 6  dataSetBuilder.AddPoint(1.4, 0.5, 0.0);
 7  dataSetBuilder.AddPoint(1.8, 0.3, 0.0);
 8  dataSetBuilder.AddPoint(0.4, 1.0, 0.0);
 9  dataSetBuilder.AddPoint(1.0, 1.2, 0.0);
10  dataSetBuilder.AddPoint(1.5, 0.9, 0.0);
11
12  dataSetBuilder.AddCell(vtkm::CELL_SHAPE_TRIANGLE);
13  dataSetBuilder.AddCellPoint(0);
14  dataSetBuilder.AddCellPoint(2);
15  dataSetBuilder.AddCellPoint(1);
16
17  dataSetBuilder.AddCell(vtkm::CELL_SHAPE_QUAD);
18  dataSetBuilder.AddCellPoint(0);
19  dataSetBuilder.AddCellPoint(4);
20  dataSetBuilder.AddCellPoint(3);
21  dataSetBuilder.AddCellPoint(2);
22
23  dataSetBuilder.AddCell(vtkm::CELL_SHAPE_TRIANGLE);
24  dataSetBuilder.AddCellPoint(1);
25  dataSetBuilder.AddCellPoint(2);
26  dataSetBuilder.AddCellPoint(5);
27
28  dataSetBuilder.AddCell(vtkm::CELL_SHAPE_POLYGON);
29  dataSetBuilder.AddCellPoint(2);
30  dataSetBuilder.AddCellPoint(3);
31  dataSetBuilder.AddCellPoint(7);
32  dataSetBuilder.AddCellPoint(6);
33  dataSetBuilder.AddCellPoint(5);
34
35  dataSetBuilder.AddCell(vtkm::CELL_SHAPE_TRIANGLE);
36  dataSetBuilder.AddCellPoint(3);
37  dataSetBuilder.AddCellPoint(4);
38  dataSetBuilder.AddCellPoint(7);
39
40  vtkm::cont::DataSet dataSet = dataSetBuilder.Create();

2.4.1.4. Add Fields

In addition to creating the geometric structure of a data set, it is usually important to add fields to the data. Fields describe numerical data associated with the topological elements in a cell. They often represent a physical quantity (such as temperature, mass, or volume fraction) but can also represent other information (such as indices or classifications).

The easiest way to define fields in a data set is to use the vtkm::cont::DataSet::AddPointField() and vtkm::cont::DataSet::AddCellField() methods. Each of these methods take a requisite field name and the array with with field data.

Both vtkm::cont::DataSet::AddPointField() and vtkm::cont::DataSet::AddCellField() are overloaded to accept arrays of data in different structures. Field arrays can be passed as standard C arrays or as std::vector’s, in which case the data are copied. Field arrays can also be passed in a ArrayHandle (introduced later in this book), in which case the data are not copied.

inline void vtkm::cont::DataSet::AddPointField(const std::string &fieldName, const vtkm::cont::UnknownArrayHandle &field)

Adds a point field of a given name to the DataSet.

Note that the indexing of fields is not the same as the order in which they are added, and that adding a field can arbitrarily reorder the integer indexing of all the fields. To retrieve a specific field, retrieve the field by name, not by integer index.

template<typename T>
inline void vtkm::cont::DataSet::AddPointField(const std::string &fieldName, const std::vector<T> &field)

Adds a point field of a given name to the DataSet.

Note that the indexing of fields is not the same as the order in which they are added, and that adding a field can arbitrarily reorder the integer indexing of all the fields. To retrieve a specific field, retrieve the field by name, not by integer index.

template<typename T>
inline void vtkm::cont::DataSet::AddPointField(const std::string &fieldName, const T *field, const vtkm::Id &n)

Adds a point field of a given name to the DataSet.

Note that the indexing of fields is not the same as the order in which they are added, and that adding a field can arbitrarily reorder the integer indexing of all the fields. To retrieve a specific field, retrieve the field by name, not by integer index.

inline void vtkm::cont::DataSet::AddCellField(const std::string &fieldName, const vtkm::cont::UnknownArrayHandle &field)

Adds a cell field of a given name to the DataSet.

Note that the indexing of fields is not the same as the order in which they are added, and that adding a field can arbitrarily reorder the integer indexing of all the fields. To retrieve a specific field, retrieve the field by name, not by integer index.

template<typename T>
inline void vtkm::cont::DataSet::AddCellField(const std::string &fieldName, const std::vector<T> &field)

Adds a cell field of a given name to the DataSet.

Note that the indexing of fields is not the same as the order in which they are added, and that adding a field can arbitrarily reorder the integer indexing of all the fields. To retrieve a specific field, retrieve the field by name, not by integer index.

template<typename T>
inline void vtkm::cont::DataSet::AddCellField(const std::string &fieldName, const T *field, const vtkm::Id &n)

Adds a cell field of a given name to the DataSet.

Note that the indexing of fields is not the same as the order in which they are added, and that adding a field can arbitrarily reorder the integer indexing of all the fields. To retrieve a specific field, retrieve the field by name, not by integer index.

The following (somewhat contrived) example defines fields for a uniform grid that identify which points and cells are on the boundary of the mesh.

Example 2.8 Adding fields to a vtkm::cont::DataSet.
 1  // Make a simple structured data set.
 2  const vtkm::Id3 pointDimensions(20, 20, 10);
 3  const vtkm::Id3 cellDimensions = pointDimensions - vtkm::Id3(1, 1, 1);
 4  vtkm::cont::DataSetBuilderUniform dataSetBuilder;
 5  vtkm::cont::DataSet dataSet = dataSetBuilder.Create(pointDimensions);
 6
 7  // Create a field that identifies points on the boundary.
 8  std::vector<vtkm::UInt8> boundaryPoints;
 9  for (vtkm::Id zIndex = 0; zIndex < pointDimensions[2]; zIndex++)
10  {
11    for (vtkm::Id yIndex = 0; yIndex < pointDimensions[1]; yIndex++)
12    {
13      for (vtkm::Id xIndex = 0; xIndex < pointDimensions[0]; xIndex++)
14      {
15        if ((xIndex == 0) || (xIndex == pointDimensions[0] - 1) || (yIndex == 0) ||
16            (yIndex == pointDimensions[1] - 1) || (zIndex == 0) ||
17            (zIndex == pointDimensions[2] - 1))
18        {
19          boundaryPoints.push_back(1);
20        }
21        else
22        {
23          boundaryPoints.push_back(0);
24        }
25      }
26    }
27  }
28
29  dataSet.AddPointField("boundary_points", boundaryPoints);
30
31  // Create a field that identifies cells on the boundary.
32  std::vector<vtkm::UInt8> boundaryCells;
33  for (vtkm::Id zIndex = 0; zIndex < cellDimensions[2]; zIndex++)
34  {
35    for (vtkm::Id yIndex = 0; yIndex < cellDimensions[1]; yIndex++)
36    {
37      for (vtkm::Id xIndex = 0; xIndex < cellDimensions[0]; xIndex++)
38      {
39        if ((xIndex == 0) || (xIndex == cellDimensions[0] - 1) || (yIndex == 0) ||
40            (yIndex == cellDimensions[1] - 1) || (zIndex == 0) ||
41            (zIndex == cellDimensions[2] - 1))
42        {
43          boundaryCells.push_back(1);
44        }
45        else
46        {
47          boundaryCells.push_back(0);
48        }
49      }
50    }
51  }
52
53  dataSet.AddCellField("boundary_cells", boundaryCells);

2.4.2. Cell Sets

A cell set determines the topological structure of the data in a data set.

class CellSet

Defines the topological structure of the data in a DataSet.

Fundamentally, any cell set is a collection of cells, which typically (but not always) represent some region in space.

Subclassed by vtkm::cont::CellSetExplicit< vtkm::cont::ArrayHandleConstant< vtkm::UInt8 >::StorageTag, ::vtkm::cont::StorageTagBasic, vtkm::cont::ArrayHandleCounting< vtkm::Id >::StorageTag >, vtkm::cont::CellSetExplicit< ShapesStorageTag, ConnectivityStorageTag, OffsetsStorageTag >, vtkm::cont::CellSetExtrude, vtkm::cont::CellSetPermutation< OriginalCellSetType_, PermutationArrayHandleType_ >, vtkm::cont::CellSetStructured< DIMENSION >

Public Functions

virtual vtkm::Id GetNumberOfCells() const = 0

Get the number of cells in the topology.

virtual vtkm::Id GetNumberOfPoints() const = 0

Get the number of points in the topology.

virtual vtkm::UInt8 GetCellShape(vtkm::Id id) const = 0

Get the shell shape of a particular cell.

virtual vtkm::IdComponent GetNumberOfPointsInCell(vtkm::Id id) const = 0

Get the number of points incident to a particular cell.

virtual void GetCellPointIds(vtkm::Id id, vtkm::Id *ptids) const = 0

Get a list of points incident to a particular cell.

virtual std::shared_ptr<CellSet> NewInstance() const = 0

Return a new CellSet that is the same derived class.

virtual void DeepCopy(const CellSet *src) = 0

Copy the provided CellSet into this object.

The provided CellSet must be the same type as this one.

virtual void PrintSummary(std::ostream&) const = 0

Print a summary of this cell set.

virtual void ReleaseResourcesExecution() = 0

Remove the CellSet from any devices.

Any memory used on a device to store this object will be deleted. However, the data will still remain on the host.

3D cells are made up of points, edges, and faces. (2D cells have only points and edges, and 1D cells have only points.) Figure 2.3 shows the relationship between a cell’s shape and these topological elements. The arrangement of these points, edges, and faces is defined by the shape of the cell, which prescribes a specific ordering of each. The basic cell shapes provided by VTK‑m are discussed in detail in Chapter 4.7 (Working with Cells).

_images/CellConstituents.png

Figure 2.3 The relationship between a cell shape and its topological elements (points, edges, and faces).

There are multiple ways to express the connections of a cell set, each with different benefits and restrictions. These different cell set types are managed by different cell set classes in VTK‑m. All VTK‑m cell set classes inherit from vtkm::cont::CellSet. The two basic types of cell sets are structured and explicit, and there are several variations of these types.

2.4.2.1. Structured Cell Sets

template<vtkm::IdComponent DIMENSION>
class CellSetStructured : public vtkm::cont::CellSet

Defines a 1-, 2-, or 3-dimensional structured grid of points.

The structured cells form lines, quadrilaterals, or hexahedra for 1-, 2-, or 3-dimensions, respectively, to connect th epoints. The topology is specified by simply providing the dimensions, which is the number of points in the i, j, and k directions of the grid of points.

Public Functions

inline virtual vtkm::Id GetNumberOfCells() const override

Get the number of cells in the topology.

inline virtual vtkm::Id GetNumberOfPoints() const override

Get the number of points in the topology.

inline virtual void ReleaseResourcesExecution() override

Remove the CellSet from any devices.

Any memory used on a device to store this object will be deleted. However, the data will still remain on the host.

inline void SetPointDimensions(SchedulingRangeType dimensions)

Set the dimensions of the structured array of points.

inline SchedulingRangeType GetPointDimensions() const

Get the dimensions of the points.

inline virtual vtkm::IdComponent GetNumberOfPointsInCell(vtkm::Id = 0) const override

Get the number of points incident to a particular cell.

inline virtual vtkm::UInt8 GetCellShape(vtkm::Id = 0) const override

Get the shell shape of a particular cell.

inline virtual void GetCellPointIds(vtkm::Id id, vtkm::Id *ptids) const override

Get a list of points incident to a particular cell.

inline virtual std::shared_ptr<CellSet> NewInstance() const override

Return a new CellSet that is the same derived class.

inline virtual void DeepCopy(const CellSet *src) override

Copy the provided CellSet into this object.

The provided CellSet must be the same type as this one.

inline virtual void PrintSummary(std::ostream &out) const override

Print a summary of this cell set.

The number of points in a vtkm::cont::CellSetStructured is implicitly \(i \times j \times k\) and the number of cells is implicitly \((i-1) \times (j-1) \times (k-1)\) (for 3D grids). Figure 2.4 demonstrates this arrangement.

_images/StructuredCellSet.png

Figure 2.4 The arrangement of points and cells in a 3D structured grid.

The big advantage of using vtkm::cont::CellSetStructured to define a cell set is that it is very space efficient because the entire topology can be defined by the three integers specifying the dimensions. Also, algorithms can be optimized for vtkm::cont::CellSetStructured’s regular nature. However, vtkm::cont::CellSetStructured’s strictly regular grid also limits its applicability. A structured cell set can only be a dense grid of lines, quadrilaterals, or hexahedra. It cannot represent irregular data well.

Many data models in other software packages, such as the one for VTK, make a distinction between uniform, rectilinear, and curvilinear grids. VTK‑m’s cell sets do not. All three of these grid types are represented by vtkm::cont::CellSetStructured. This is because in a VTK‑m data set the cell set and the coordinate system are defined independently and used interchangeably. A structured cell set with uniform point coordinates makes a uniform grid. A structured cell set with point coordinates defined irregularly along coordinate axes makes a rectilinear grid. And a structured cell set with arbitrary point coordinates makes a curvilinear grid. The point coordinates are defined by the data set’s coordinate system, which is discussed in Section 2.4.4 (Coordinate Systems).

2.4.2.2. Explicit Cell Sets

template<typename ShapesStorageTag = ::vtkm::cont::StorageTagBasic, typename ConnectivityStorageTag = ::vtkm::cont::StorageTagBasic, typename OffsetsStorageTag = ::vtkm::cont::StorageTagBasic>
class CellSetExplicit : public vtkm::cont::CellSet

Defines an irregular collection of cells.

The cells can be of different types and connected in arbitrary ways. This is done by explicitly providing for each cell a sequence of points that defines the cell.

Public Functions

virtual vtkm::Id GetNumberOfCells() const override

Get the number of cells in the topology.

virtual vtkm::Id GetNumberOfPoints() const override

Get the number of points in the topology.

virtual void PrintSummary(std::ostream &out) const override

Print a summary of this cell set.

virtual void ReleaseResourcesExecution() override

Remove the CellSet from any devices.

Any memory used on a device to store this object will be deleted. However, the data will still remain on the host.

virtual std::shared_ptr<CellSet> NewInstance() const override

Return a new CellSet that is the same derived class.

virtual void DeepCopy(const CellSet *src) override

Copy the provided CellSet into this object.

The provided CellSet must be the same type as this one.

virtual vtkm::IdComponent GetNumberOfPointsInCell(vtkm::Id cellid) const override

Get the number of points incident to a particular cell.

virtual void GetCellPointIds(vtkm::Id id, vtkm::Id *ptids) const override

Get a list of points incident to a particular cell.

virtual vtkm::UInt8 GetCellShape(vtkm::Id cellid) const override

Get the shell shape of a particular cell.

void PrepareToAddCells(vtkm::Id numCells, vtkm::Id connectivityMaxLen)

Start adding cells one at a time.

After this method is called, AddCell is called repeatedly to add each cell. Once all cells are added, call CompleteAddingCells.

template<typename IdVecType>
void AddCell(vtkm::UInt8 cellType, vtkm::IdComponent numVertices, const IdVecType &ids)

Add a cell.

This can only be called after AddCell.

void CompleteAddingCells(vtkm::Id numPoints)

Finish adding cells one at a time.

void Fill(vtkm::Id numPoints, const vtkm::cont::ArrayHandle<vtkm::UInt8, ShapesStorageTag> &cellTypes, const vtkm::cont::ArrayHandle<vtkm::Id, ConnectivityStorageTag> &connectivity, const vtkm::cont::ArrayHandle<vtkm::Id, OffsetsStorageTag> &offsets)

Set all the cells of the mesh.

This method can be used to fill the memory from another system without copying data.

The types of cell sets are listed in Figure 2.5.

_images/CellConnections.png

Figure 2.5 Basic Cell Shapes in a vtkm::cont::CellSetExplicit.

An explicit cell set is defined with a minimum of three arrays. The first array identifies the shape of each cell. (Identifiers for cell shapes are shown in Figure 2.5.) The second array has a sequence of point indices that make up each cell. The third array identifies an offset into the second array where the point indices for each cell is found plus an extra entry at the end set to the size of the second array. Figure 2.6 shows a simple example of an explicit cell set.

_images/ExplicitCellConnections.png

Figure 2.6 Example of cells in a vtkm::cont::CellSetExplicit and the arrays that define them.

An explicit cell set can also identify the number of indices defined for each cell by subtracting consecutive entries in the offsets array. It is often the case when creating a vtkm::cont::CellSetExplicit that you have an array containing the number of indices rather than the offsets. Such an array can be converted to an offsets array that can be used with vtkm::cont::CellSetExplicit by using the vtkm::cont::ConvertNumComponentsToOffsets() convenience function.

void vtkm::cont::ConvertNumComponentsToOffsets(const vtkm::cont::UnknownArrayHandle &numComponentsArray, vtkm::cont::ArrayHandle<vtkm::Id> &offsetsArray, vtkm::Id &componentsArraySize, vtkm::cont::DeviceAdapterId device = vtkm::cont::DeviceAdapterTagAny{})

ConvertNumComponentsToOffsets takes an array of Vec sizes (i.e.

the number of components in each Vec) and returns an array of offsets to a packed array of such Vecs. The resulting array can be used with ArrayHandleGroupVecVariable.

Note that this function is pre-compiled for some set of ArrayHandle types. If you get a warning about an inefficient conversion (or the operation fails outright), you might need to use vtkm::cont::internal::ConvertNumComponentsToOffsetsTemplate.

Parameters:
  • numComponentsArray[in] the input array that specifies the number of components in each group Vec.

  • offsetsArray[out] (optional) the output ArrayHandle, which must have a value type of vtkm::Id. If the output ArrayHandle is not given, it is returned.

  • componentsArraySize[in] (optional) a reference to a vtkm::Id and is filled with the expected size of the component values array.

  • device[in] (optional) specifies the device on which to run the conversion.

vtkm::cont::CellSetExplicit is a powerful representation for a cell set because it can represent an arbitrary collection of cells. However, because all connections must be explicitly defined, vtkm::cont::CellSetExplicit requires a significant amount of memory to represent the topology.

An important specialization of an explicit cell set is vtkm::cont::CellSetSingleType.

template<typename ConnectivityStorageTag = ::vtkm::cont::StorageTagBasic>
class CellSetSingleType : public vtkm::cont::CellSetExplicit<vtkm::cont::ArrayHandleConstant<vtkm::UInt8>::StorageTag, ::vtkm::cont::StorageTagBasic, vtkm::cont::ArrayHandleCounting<vtkm::Id>::StorageTag>

An explicit cell set with all cells of the same shape.

CellSetSingleType is an explicit cell set constrained to contain cells that all have the same shape and all have the same number of points. So, for example if you are creating a surface that you know will contain only triangles, CellSetSingleType is a good representation for these data.

Using CellSetSingleType saves memory because the array of cell shapes and the array of point counts no longer need to be stored. CellSetSingleType also allows VTK-m to skip some processing and other storage required for general explicit cell sets.

Public Functions

inline void PrepareToAddCells(vtkm::Id numCells, vtkm::Id connectivityMaxLen)

Start adding cells one at a time.

After this method is called, AddCell is called repeatedly to add each cell. Once all cells are added, call CompleteAddingCells.

template<typename IdVecType>
inline void AddCell(vtkm::UInt8 shapeId, vtkm::IdComponent numVertices, const IdVecType &ids)

Add a cell.

This can only be called after AddCell.

inline void CompleteAddingCells(vtkm::Id numPoints)

Finish adding cells one at a time.

inline void Fill(vtkm::Id numPoints, vtkm::UInt8 shapeId, vtkm::IdComponent numberOfPointsPerCell, const vtkm::cont::ArrayHandle<vtkm::Id, ConnectivityStorageTag> &connectivity)

Set all the cells of the mesh.

This method can be used to fill the memory from another system without copying data.

inline virtual vtkm::UInt8 GetCellShape(vtkm::Id) const override

Get the shell shape of a particular cell.

inline virtual std::shared_ptr<CellSet> NewInstance() const override

Return a new CellSet that is the same derived class.

inline virtual void DeepCopy(const CellSet *src) override

Copy the provided CellSet into this object.

The provided CellSet must be the same type as this one.

inline virtual void PrintSummary(std::ostream &out) const override

Print a summary of this cell set.

2.4.2.3. Cell Set Permutations

To rearrange, and possibly subsample, cells in a CellSet, use vtkm::cont::CellSetPermutation to define a new set without copying.

template<typename OriginalCellSetType_, typename PermutationArrayHandleType_ = vtkm::cont::ArrayHandle<vtkm::Id, ::vtkm::cont::StorageTagBasic>>
class CellSetPermutation : public vtkm::cont::CellSet

Rearranges the cells of one cell set to create another cell set.

This restructuring of cells is not done by copying data to a new structure. Rather, CellSetPermutation establishes a look-up from one cell structure to another. Cells are permuted on the fly while algorithms are run.

A CellSetPermutation is established by providing a mapping array that for every cell index provides the equivalent cell index in the cell set being permuted. CellSetPermutation is most often used to mask out cells in a data set so that algorithms will skip over those cells when running.

Public Functions

inline CellSetPermutation(const PermutationArrayHandleType &validCellIds, const OriginalCellSetType &cellset)

Create a CellSetPermutation.

Parameters:
  • validCellIds[in] An array that defines the permutation. If index i is value j, then the ith cell of this cell set will be the same as the jth cell in the original cellset.

  • cellset[in] The original cell set that this one is permuting.

inline const OriginalCellSetType &GetFullCellSet() const

Returns the original CellSet that this one is permuting.

inline const PermutationArrayHandleType &GetValidCellIds() const

Returns the array used to permute the cell indices.

inline virtual vtkm::Id GetNumberOfCells() const override

Get the number of cells in the topology.

inline virtual vtkm::Id GetNumberOfPoints() const override

Get the number of points in the topology.

inline virtual void ReleaseResourcesExecution() override

Remove the CellSet from any devices.

Any memory used on a device to store this object will be deleted. However, the data will still remain on the host.

inline virtual vtkm::IdComponent GetNumberOfPointsInCell(vtkm::Id cellIndex) const override

Get the number of points incident to a particular cell.

inline virtual vtkm::UInt8 GetCellShape(vtkm::Id id) const override

Get the shell shape of a particular cell.

inline virtual void GetCellPointIds(vtkm::Id id, vtkm::Id *ptids) const override

Get a list of points incident to a particular cell.

inline virtual std::shared_ptr<CellSet> NewInstance() const override

Return a new CellSet that is the same derived class.

inline virtual void DeepCopy(const CellSet *src) override

Copy the provided CellSet into this object.

The provided CellSet must be the same type as this one.

inline void Fill(const PermutationArrayHandleType &validCellIds, const OriginalCellSetType &cellset)

Set the topology.

Parameters:
  • validCellIds[in] An array that defines the permutation. If index i is value j, then the ith cell of this cell set will be the same as the jth cell in the original cellset.

  • cellset[in] The original cell set that this one is permuting.

inline virtual void PrintSummary(std::ostream &out) const override

Print a summary of this cell set.

Did You Know?

Although vtkm::cont::CellSetPermutation can mask cells, it cannot mask points. All points from the original cell set are available in the permuted cell set regardless of whether they are used.

The following example uses vtkm::cont::CellSetPermutation with a counting array to expose every tenth cell. This provides a simple way to subsample a data set.

Example 2.9 Subsampling a data set with vtkm::cont::CellSetPermutation.
 1  // Create a simple data set.
 2  vtkm::cont::DataSetBuilderUniform dataSetBuilder;
 3  vtkm::cont::DataSet originalDataSet = dataSetBuilder.Create(vtkm::Id3(33, 33, 26));
 4  vtkm::cont::CellSetStructured<3> originalCellSet;
 5  originalDataSet.GetCellSet().AsCellSet(originalCellSet);
 6
 7  // Create a permutation array for the cells. Each value in the array refers
 8  // to a cell in the original cell set. This particular array selects every
 9  // 10th cell.
10  vtkm::cont::ArrayHandleCounting<vtkm::Id> permutationArray(0, 10, 2560);
11
12  // Create a permutation of that cell set containing only every 10th cell.
13  vtkm::cont::CellSetPermutation<vtkm::cont::CellSetStructured<3>,
14                                 vtkm::cont::ArrayHandleCounting<vtkm::Id>>
15    permutedCellSet(permutationArray, originalCellSet);

2.4.2.4. Cell Set Extrude

class CellSetExtrude : public vtkm::cont::CellSet

Defines a 3-dimensional extruded mesh representation.

CellSetExtrude takes takes a mesh defined in the XZ-plane and extrudes it along the Y-axis. This plane is repeated in a series of steps and forms wedge cells between them.

The extrusion can be linear or rotational (e.g., to form a torus).

Public Functions

virtual vtkm::Id GetNumberOfCells() const override

Get the number of cells in the topology.

virtual vtkm::Id GetNumberOfPoints() const override

Get the number of points in the topology.

virtual vtkm::UInt8 GetCellShape(vtkm::Id id) const override

Get the shell shape of a particular cell.

virtual vtkm::IdComponent GetNumberOfPointsInCell(vtkm::Id id) const override

Get the number of points incident to a particular cell.

virtual void GetCellPointIds(vtkm::Id id, vtkm::Id *ptids) const override

Get a list of points incident to a particular cell.

virtual std::shared_ptr<CellSet> NewInstance() const override

Return a new CellSet that is the same derived class.

virtual void DeepCopy(const CellSet *src) override

Copy the provided CellSet into this object.

The provided CellSet must be the same type as this one.

virtual void PrintSummary(std::ostream &out) const override

Print a summary of this cell set.

virtual void ReleaseResourcesExecution() override

Remove the CellSet from any devices.

Any memory used on a device to store this object will be deleted. However, the data will still remain on the host.

_images/ExtrudedCellSet.png

Figure 2.7 An example of an extruded wedge from XZ-plane coordinates. Six wedges are extracted from three XZ-plane points.

The extruded mesh is advantageous because it is represented on-the-fly as required, so no additional memory is required. In contrast other forms of cell sets, such as vtkm::cont::CellSetExplicit, need to be explicitly constructed by replicating the vertices and cells. Figure 2.7 shows an example of six wedges extruded from three 2-dimensional coordinates.

2.4.2.5. Unknown Cell Sets

Each of the aforementioned cell set types are represented by a different class. A vtkm::cont::DataSet object must hold one of these cell set objects that represent the cell structure. The actual object used is not determined until run time.

The vtkm::cont::DataSet object manages the cell set object with vtkm::cont::UnknownCellSet. When you call vtkm::cont::DataSet::GetCellSet(), it returns a vtkm::cont::UnknownCellSet.

The vtkm::cont::UnknownCellSet object provides mechanisms to query the cell set, identify its type, and cast it to one of the concrete CellSet types. See Chapter ref{chap:UnknownCellSet} for details on working with vtkm::cont::UnknownCellSet.

2.4.3. Fields

A field on a data set provides a value on every point in space on the mesh. Fields are often used to describe physical properties such as pressure, temperature, mass, velocity, and much more. Fields are represented in a VTK‑m data set as an array where each value is associated with a particular element type of a mesh (such as points or cells). This association of field values to mesh elements and the structure of the cell set determines how the field is interpolated throughout the space of the mesh.

Fields are manged by the vtkm::cont::Field class.

class Field

A Field encapsulates an array on some piece of the mesh, such as the points, a cell set, a point logical dimension, or the whole mesh.

Subclassed by vtkm::cont::CoordinateSystem

Fields are identified by a simple name string.

inline const std::string &vtkm::cont::Field::GetName() const

The vtkm::cont::Field object internally holds a reference to an array in a type-agnostic way. Filters and other VTK‑m units will determine the type of the array and pull it out of the vtkm::cont::Field.

const vtkm::cont::UnknownArrayHandle &vtkm::cont::Field::GetData() const

The field data is associated with a particular type of element of a mesh such as points, cells, or the whole mesh.

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

Associations are identified by the vtkm::cont::Field::Association enumeration.

enum class vtkm::cont::Field::Association

Identifies what elements of a data set a field is associated with.

The Association enum is used by vtkm::cont::Field to specify on what topological elements each item in the field is associated with.

Values:

enumerator Any

Any field regardless of the association.

This is used when choosing a vtkm::cont::Field that could be of any association. It is often used as the default if no association is given.

enumerator WholeDataSet

A “global” field that applies to the entirety of a vtkm::cont::DataSet.

Fields of this association often contain summary or annotation information. An example of a whole data set field could be the region that the mesh covers.

enumerator Points

A field that applies to points.

There is a separate field value attached to each point. Point fields usually represent samples of continuous data that can be reinterpolated through cells. Physical properties such as temperature, pressure, density, velocity, etc. are usually best represented in point fields. Data that deals with the points of the topology, such as displacement vectors, are also appropriate for point data.

enumerator Cells

A field that applies to cells.

There is a separate field value attached to each cell in a cell set. Cell fields usually represent values from an integration over the finite cells of the mesh. Integrated values like mass or volume are best represented in cell fields. Statistics about each cell like strain or cell quality are also appropriate for cell data.

enumerator Partitions

A field that applies to partitions.

This type of field is attached to a vtkm::cont::PartitionedDataSet. There is a separate field value attached to each partition. Identification or information about the arrangement of partitions such as hierarchy levels are usually best represented in partition fields.

enumerator Global

A field that applies to all partitions.

This type of field is attached to a vtkm::cont::PartitionedDataSet. It contains values that are “global” across all partitions and data therin.

The vtkm::cont::Field class also has several convenience methods for querying the association.

inline bool vtkm::cont::Field::IsPointField() const
inline bool vtkm::cont::Field::IsCellField() const
inline bool vtkm::cont::Field::IsWholeDataSetField() const
inline bool vtkm::cont::Field::IsPartitionsField() const
inline bool vtkm::cont::Field::IsGlobalField() const

vtkm::cont::Field has a convenience method named vtkm::cont::Field::GetRange() that finds the range of values stored in the field array.

const vtkm::cont::ArrayHandle<vtkm::Range> &vtkm::cont::Field::GetRange() const

Returns the range of each component in the field array.

The ranges of each component are returned in an ArrayHandle containing vtkm::Range values. So, for example, calling GetRange on a scalar field will return an ArrayHandle with exactly 1 entry in it. Calling GetRange on a field of 3D vectors will return an ArrayHandle with exactly 3 entries corresponding to each of the components in the range.

Details on how to get data from a vtkm::cont::ArrayHandle them is given in Chapter ref{chap:AccessingAllocatingArrays}.

2.4.4. Coordinate Systems

A coordinate system determines the location of a mesh’s elements in space. The spatial location is described by providing a 3D vector at each point that gives the coordinates there. The point coordinates can then be interpolated throughout the mesh.

class CoordinateSystem : public vtkm::cont::Field

Manages a coordinate system for a DataSet.

A coordinate system is really a field with a special meaning, so CoordinateSystem class inherits from the Field class. CoordinateSystem constrains the field to be associated with points and typically has 3D floating point vectors for values.

In addition to all the methods provided by the vtkm::cont::Field superclass, the vtkm::cont::CoordinateSystem also provides a vtkm::cont::CoordinateSystem::GetBounds() convenience method that returns a vtkm::Bounds object giving the spatial bounds of the coordinate system.

inline vtkm::Bounds vtkm::cont::CoordinateSystem::GetBounds() const

It is typical for a vtkm::cont::DataSet to have one coordinate system defined, but it is possible to define multiple coordinate systems. This is helpful when there are multiple ways to express coordinates. For example, positions in geographic may be expressed as Cartesian coordinates or as latitude-longitude coordinates. Both are valid and useful in different ways.

It is also valid to have a vtkm::cont::DataSet with no coordinate system. This is useful when the structure is not rooted in physical space. For example, if the cell set is representing a graph structure, there might not be any physical space that has meaning for the graph.

2.4.5. Partitioned Data Sets

class PartitionedDataSet

Comprises a set of vtkm::cont::DataSet objects.

Iterators

PartitionedDataSet provides an iterator interface that allows you to iterate over the contained partitions using the for (auto ds : pds) syntax.

Public Functions

PartitionedDataSet(const vtkm::cont::DataSet &ds)

Create a new PartitionedDataSet containng a single DataSet ds.

explicit PartitionedDataSet(const std::vector<vtkm::cont::DataSet> &partitions)

Create a new PartitionedDataSet with a DataSet vector partitions.

explicit PartitionedDataSet(vtkm::Id size)

Create a new PartitionedDataSet with the capacity set to be size.

vtkm::cont::Field GetFieldFromPartition(const std::string &field_name, int partition_index) const

Get the field field_name from partition partition_index.

vtkm::Id GetNumberOfPartitions() const

Get number of DataSet objects stored in this PartitionedDataSet.

vtkm::Id GetGlobalNumberOfPartitions() const

Get number of partations across all MPI ranks.

Warning

This method requires global communication (MPI_Allreduce) if MPI is enabled.

const vtkm::cont::DataSet &GetPartition(vtkm::Id partId) const

Get the DataSet partId.

const std::vector<vtkm::cont::DataSet> &GetPartitions() const

Get an STL vector of all DataSet objects stored in PartitionedDataSet.

void AppendPartition(const vtkm::cont::DataSet &ds)

Add DataSet ds to the end of the list of partitions.

void InsertPartition(vtkm::Id index, const vtkm::cont::DataSet &ds)

Add DataSet ds to position index of the contained DataSet vector.

All partitions at or after this location are pushed back.

void ReplacePartition(vtkm::Id index, const vtkm::cont::DataSet &ds)

Replace the index positioned element of the contained DataSet vector with ds.

void AppendPartitions(const std::vector<vtkm::cont::DataSet> &partitions)

Append the DataSet vector partitions to the end of list of partitions.

This list can be provided as a std::vector, or it can be an initializer list (declared in { } curly braces).

inline vtkm::IdComponent GetNumberOfFields() const

Methods to Add and Get fields on a PartitionedDataSet.

inline void AddField(const Field &field)

Adds a field that is applied to the meta-partition structure.

The field must have a partition that applies across all partitions.

template<typename T, typename Storage>
inline void AddGlobalField(const std::string &fieldName, const vtkm::cont::ArrayHandle<T, Storage> &field)

Add a field with a global association.

template<typename T, typename Storage>
inline void AddPartitionsField(const std::string &fieldName, const vtkm::cont::ArrayHandle<T, Storage> &field)

Add a field where each entry is associated with a whole partition.

inline vtkm::cont::Field &GetField(const std::string &name, vtkm::cont::Field::Association assoc = vtkm::cont::Field::Association::Any)

Get a field associated with the partitioned data structure.

The field is selected by name and, optionally, the association.

inline const vtkm::cont::Field &GetGlobalField(const std::string &name) const

Get a global field.

inline const vtkm::cont::Field &GetPartitionsField(const std::string &name) const

Get a field associated with the partitions.

inline bool HasField(const std::string &name, vtkm::cont::Field::Association assoc = vtkm::cont::Field::Association::Any) const

Query whether the partitioned data set has the named field.

inline bool HasGlobalField(const std::string &name) const

Query whether the partitioned data set has the named global field.

inline bool HasPartitionsField(const std::string &name) const

Query whether the partitioned data set has the named partition field.

void CopyPartitions(const vtkm::cont::PartitionedDataSet &source)

Copies the partitions from the source. The fields on the PartitionedDataSet are not copied.

The following example creates a vtkm::cont::PartitionedDataSet containing two uniform grid data sets.

Example 2.10 Creating a vtkm::cont::PartitionedDataSet.
1  // Create two uniform data sets
2  vtkm::cont::DataSetBuilderUniform dataSetBuilder;
3
4  vtkm::cont::DataSet dataSet1 = dataSetBuilder.Create(vtkm::Id3(10, 10, 10));
5  vtkm::cont::DataSet dataSet2 = dataSetBuilder.Create(vtkm::Id3(30, 30, 30));
6
7  // Add the datasets to a multi block
8  vtkm::cont::PartitionedDataSet partitionedData;
9  partitionedData.AppendPartitions({ dataSet1, dataSet2 });

It is always possible to retrieve the independent blocks in a vtkm::cont::PartitionedDataSet, from which you can iterate and get information about the data. However, VTK‑m provides several helper functions to collect metadata information about the collection as a whole.

vtkm::Bounds vtkm::cont::BoundsCompute(const vtkm::cont::DataSet &dataset, vtkm::Id coordinate_system_index = 0)

Functions to compute bounds for a single dataSset or partition dataset.

These are utility functions that compute bounds for a single dataset or partitioned dataset. When VTK-m is operating in an distributed environment, these are bounds on the local process. To get global bounds across all ranks, use vtkm::cont::BoundsGlobalCompute instead.

Note that if the provided CoordinateSystem does not exists, empty bounds are returned. Likewise, for PartitionedDataSet, partitions without the chosen CoordinateSystem are skipped.

vtkm::Bounds vtkm::cont::BoundsCompute(const vtkm::cont::PartitionedDataSet &pds, vtkm::Id coordinate_system_index = 0)
vtkm::Bounds vtkm::cont::BoundsCompute(const vtkm::cont::DataSet &dataset, const std::string &coordinate_system_name)
vtkm::Bounds vtkm::cont::BoundsCompute(const vtkm::cont::PartitionedDataSet &pds, const std::string &coordinate_system_name)
vtkm::Bounds vtkm::cont::BoundsGlobalCompute(const vtkm::cont::DataSet &dataset, vtkm::Id coordinate_system_index = 0)

Functions to compute bounds for a single dataset or partitioned dataset globally.

These are utility functions that compute bounds for a single dataset or partitioned dataset globally i.e. across all ranks when operating in a distributed environment. When VTK-m not operating in an distributed environment, these behave same as vtkm::cont::BoundsCompute.

Note that if the provided CoordinateSystem does not exists, empty bounds are returned. Likewise, for PartitionedDataSet, partitions without the chosen CoordinateSystem are skipped.

vtkm::Bounds vtkm::cont::BoundsGlobalCompute(const vtkm::cont::PartitionedDataSet &pds, vtkm::Id coordinate_system_index = 0)
vtkm::Bounds vtkm::cont::BoundsGlobalCompute(const vtkm::cont::DataSet &dataset, const std::string &coordinate_system_name)
vtkm::Bounds vtkm::cont::BoundsGlobalCompute(const vtkm::cont::PartitionedDataSet &pds, const std::string &coordinate_system_name)
vtkm::cont::ArrayHandle<vtkm::Range> vtkm::cont::FieldRangeCompute(const vtkm::cont::DataSet &dataset, const std::string &name, vtkm::cont::Field::Association assoc = vtkm::cont::Field::Association::Any)

Compute ranges for fields in a DataSet or PartitionedDataSet.

These methods to compute ranges for fields in a single dataset or a partitioned dataset. When using VTK-m in a hybrid-parallel environment with distributed processing, this class uses ranges for locally available data alone. Use FieldRangeGlobalCompute to compute ranges globally across all ranks even in distributed mode. Returns the range for a field from a dataset. If the field is not present, an empty ArrayHandle will be returned.

vtkm::cont::ArrayHandle<vtkm::Range> vtkm::cont::FieldRangeCompute(const vtkm::cont::PartitionedDataSet &pds, const std::string &name, vtkm::cont::Field::Association assoc = vtkm::cont::Field::Association::Any)

Returns the range for a field from a PartitionedDataSet.

If the field is not present on any of the partitions, an empty ArrayHandle will be returned. If the field is present on some partitions, but not all, those partitions without the field are skipped.

The returned array handle will have as many values as the maximum number of components for the selected field across all partitions.

vtkm::cont::ArrayHandle<vtkm::Range> vtkm::cont::FieldRangeGlobalCompute(const vtkm::cont::DataSet &dataset, const std::string &name, vtkm::cont::Field::Association assoc = vtkm::cont::Field::Association::Any)

utility functions to compute global ranges for dataset fields.

These functions compute global ranges for fields in a single DataSet or a PartitionedDataSet. In non-distributed environments, this is exactly same as FieldRangeCompute. In distributed environments, however, the range is computed locally on each rank and then a reduce-all collective is performed to reduces the ranges on all ranks. Returns the range for a field from a dataset. If the field is not present, an empty ArrayHandle will be returned.

vtkm::cont::ArrayHandle<vtkm::Range> vtkm::cont::FieldRangeGlobalCompute(const vtkm::cont::PartitionedDataSet &pds, const std::string &name, vtkm::cont::Field::Association assoc = vtkm::cont::Field::Association::Any)

Returns the range for a field from a PartitionedDataSet.

If the field is not present on any of the partitions, an empty ArrayHandle will be returned. If the field is present on some partitions, but not all, those partitions without the field are skipped.

The returned array handle will have as many values as the maximum number of components for the selected field across all partitions.

The following example illustrates a spatial bounds query and a field range query on a vtkm::cont::PartitionedDataSet.

Example 2.11 Queries on a vtkm::cont::PartitionedDataSet.
1  // Get the bounds of a multi-block data set
2  vtkm::Bounds bounds = vtkm::cont::BoundsCompute(partitionedData);
3
4  // Get the overall min/max of a field named "cellvar"
5  vtkm::cont::ArrayHandle<vtkm::Range> cellvarRanges =
6    vtkm::cont::FieldRangeCompute(partitionedData, "cellvar");
7
8  // Assuming the "cellvar" field has scalar values, then cellvarRanges has one entry
9  vtkm::Range cellvarRange = cellvarRanges.ReadPortal().Get(0);

Did You Know?

The aforementioned functions for querying a vtkm::cont::PartitionedDataSet object also work on vtkm::cont::DataSet objects. This is particularly useful with the vtkm::cont::BoundsGlobalCompute() and vtkm::cont::FieldRangeGlobalCompute() functions to manage distributed parallel objects.

Filters can be executed on vtkm::cont::PartitionedDataSet objects in a similar way they are executed on vtkm::cont::DataSet objects. In both cases, the vtkm::cont::Filter::Execute() method is called on the filter giving data object as an argument.

Example 2.12 Applying a filter to multi block data.
1  vtkm::filter::field_conversion::CellAverage cellAverage;
2  cellAverage.SetActiveField("pointvar", vtkm::cont::Field::Association::Points);
3
4  vtkm::cont::PartitionedDataSet results = cellAverage.Execute(partitionedData);