VTK  9.3.20240420
Namespaces | Classes | Typedefs | Enumerations | Functions
vtk Namespace Reference

Specialization of tuple ranges and iterators for vtkAOSDataArrayTemplate. More...

Namespaces

namespace  detail
 
namespace  hypertreegrid
 
namespace  literals
 

Classes

class  CompositeDataSetNodeReference
 A reference proxy into a vtkCompositeDataSet, obtained by dereferencing an iterator from the vtk::Range(vtkCompositeDataSet*) overloads. More...
 
class  HasSuperclass
 Determine whether the provided class (VTKObjectType ) has a parent class. More...
 
struct  ParentClasses
 Invoke a functor on the named type and each of its parent types. More...
 
struct  ParentClasses< VTKObjectType, false >
 
struct  ParentClasses< VTKObjectType, true >
 

Typedefs

using ComponentIdType = int
 
using TupleIdType = vtkIdType
 
using ValueIdType = vtkIdType
 
template<typename ArrayType , typename ForceValueTypeForVtkDataArray = double, typename = detail::EnableIfVtkDataArray<ArrayType>>
using GetAPIType = typename detail::GetAPITypeImpl< ArrayType, ForceValueTypeForVtkDataArray >::APIType
 
template<typename ArrayType >
using IsAOSDataArray = std::integral_constant< bool, detail::IsAOSDataArrayImpl< ArrayType >::value >
 

Enumerations

enum class  CompositeDataSetOptions : unsigned int { None = 0 , SkipEmptyNodes = 1 << 1 }
 
enum class  DataObjectTreeOptions : unsigned int { None = 0 , SkipEmptyNodes = 1 << 1 , VisitOnlyLeaves = 1 << 2 , TraverseSubTree = 1 << 3 }
 

Functions

vtkSmartPointer< vtkCompositeArray< T > > ConcatenateDataArrays (const std::vector< vtkDataArray * > &arrays)
 
template<ComponentIdType TupleSize = detail::DynamicTupleSize, typename ArrayTypePtr = vtkDataArray*>
VTK_ITER_INLINE auto DataArrayTupleRange (const ArrayTypePtr &array, TupleIdType start=-1, TupleIdType end=-1) -> typename detail::SelectTupleRange< ArrayTypePtr, TupleSize >::type
 Generate an stl and for-range compatible range of tuple iterators from a vtkDataArray.
 
template<ComponentIdType TupleSize = detail::DynamicTupleSize, typename ForceValueTypeForVtkDataArray = double, typename ArrayTypePtr = vtkDataArray*>
VTK_ITER_INLINE auto DataArrayValueRange (const ArrayTypePtr &array, ValueIdType start=-1, ValueIdType end=-1) -> typename detail::SelectValueRange< ArrayTypePtr, TupleSize, ForceValueTypeForVtkDataArray >::type
 Generate an stl and for-range compatible range of flat AOS iterators from a vtkDataArray.
 
template<typename IterablePtr , typename... Options>
auto Range (IterablePtr iterable, Options &&... opts) -> typename detail::IterableTraits< typename detail::StripPointers< IterablePtr >::type >::RangeType
 Generate an iterable STL proxy object for a VTK container.
 
template<typename T >
vtkSmartPointer< T > MakeSmartPointer (T *obj)
 Construct a vtkSmartPointer<T> containing obj.
 
template<typename T >
vtkSmartPointer< T > TakeSmartPointer (T *obj)
 Construct a vtkSmartPointer<T> containing obj.
 
template<typename ValueType >
vtkSmartPointer< vtkImplicitArray< vtkStructuredPointBackend< ValueType > > > CreateStructuredPointArray (vtkDataArray *xCoords, vtkDataArray *yCoords, vtkDataArray *zCoords, int extent[6], int dataDescription, double dirMatrix[9])
 Create an implicit point array from the given coordinates and direction matrix which is optional.
 
template<typename ObjectType >
std::string TypeName ()
 Return the demangled type-name of the provided ObjectType.
 
template<typename ObjectType >
vtkStringToken TypeToken ()
 Return a string token holding a hash of the demangled type-name of the provided ObjectType.
 
template<typename T >
vtkWeakPointer< T > TakeWeakPointer (T *obj)
 Construct a vtkWeakPointer<T> containing obj.
 
template<typename VTKObjectType , typename Container >
void Inherits (Container &container)
 Populate the container with the name of this class and its ancestors.
 
template<typename VTKObjectType , typename StopAtType , typename Container >
void Inherits (Container &container)
 Populate the container with the name of this class and its ancestors.
 

Detailed Description

Specialization of tuple ranges and iterators for vtkAOSDataArrayTemplate.

Compare smart pointer values.

This file contains the traits for the implicit array mechanism in VTK.

Generic implementation of value ranges and iterators, suitable for vtkDataArray and all subclasses.

Specialization of value ranges and iterators for vtkAOSDataArrayTemplate.

Generic implementation of tuple ranges and iterators, suitable for vtkDataArray and all subclasses.

These traits are very much an internal to vtkImplicitArrays and normal developers looking to develop a new vtkImplicitArray should (ideally) not have to open this file.

In order to ensure that template parameters passed to the vtkImplicitArray share a common interface without having to subclass all of them from the same abstract class, we have decided to use a trait mechanism to statically dispatch the functionalities of types passed as template parameters to the array.

There is 1 mandatory traits that a template type to vtkImplicitArray must implement:

Potential improvements to implicit arrays which would allow for write access would include the following 2 optional traits:

All the traits defining the behavior of the implicit "function" or "backend" to the vtkImplicitArray should be composited into the implicit_array_traits

Typedef Documentation

◆ ComponentIdType

using vtk::ComponentIdType = typedef int

Definition at line 60 of file vtkDataArrayMeta.h.

◆ TupleIdType

using vtk::TupleIdType = typedef vtkIdType

Definition at line 61 of file vtkDataArrayMeta.h.

◆ ValueIdType

using vtk::ValueIdType = typedef vtkIdType

Definition at line 62 of file vtkDataArrayMeta.h.

◆ GetAPIType

template<typename ArrayType , typename ForceValueTypeForVtkDataArray = double, typename = detail::EnableIfVtkDataArray<ArrayType>>
using vtk::GetAPIType = typedef typename detail::GetAPITypeImpl<ArrayType, ForceValueTypeForVtkDataArray>::APIType

Definition at line 186 of file vtkDataArrayMeta.h.

◆ IsAOSDataArray

template<typename ArrayType >
using vtk::IsAOSDataArray = typedef std::integral_constant<bool, detail::IsAOSDataArrayImpl<ArrayType>::value>

Definition at line 210 of file vtkDataArrayMeta.h.

Enumeration Type Documentation

◆ CompositeDataSetOptions

enum class vtk::CompositeDataSetOptions : unsigned int
strong
Enumerator
None 
SkipEmptyNodes 

Definition at line 21 of file vtkCompositeDataSetRange.h.

◆ DataObjectTreeOptions

enum class vtk::DataObjectTreeOptions : unsigned int
strong
Enumerator
None 
SkipEmptyNodes 
VisitOnlyLeaves 
TraverseSubTree 

Definition at line 21 of file vtkDataObjectTreeRange.h.

Function Documentation

◆ ConcatenateDataArrays()

vtkSmartPointer< vtkCompositeArray< T > > vtk::ConcatenateDataArrays ( const std::vector< vtkDataArray * > &  arrays)

◆ DataArrayTupleRange()

template<ComponentIdType TupleSize = detail::DynamicTupleSize, typename ArrayTypePtr = vtkDataArray*>
VTK_ITER_INLINE auto vtk::DataArrayTupleRange ( const ArrayTypePtr &  array,
TupleIdType  start = -1,
TupleIdType  end = -1 
) -> typename detail::SelectTupleRange<ArrayTypePtr, TupleSize>::type

Generate an stl and for-range compatible range of tuple iterators from a vtkDataArray.

This function returns a TupleRange object that is compatible with C++11 for-range syntax. As an example usage, consider a function that takes some instance of vtkDataArray (or a subclass) and prints the magnitude of each tuple:

template <typename ArrayType>
void PrintMagnitudes(ArrayType *array)
{
for (const auto tuple : vtk::DataArrayTupleRange(array))
{
double mag = 0.;
for (const T comp : tuple)
{
mag += static_cast<double>(comp) * static_cast<double>(comp);
}
mag = std::sqrt(mag);
std::cerr << mag < "\n";
}
}
Specialization of tuple ranges and iterators for vtkAOSDataArrayTemplate.
VTK_ITER_INLINE auto DataArrayTupleRange(const ArrayTypePtr &array, TupleIdType start=-1, TupleIdType end=-1) -> typename detail::SelectTupleRange< ArrayTypePtr, TupleSize >::type
Generate an stl and for-range compatible range of tuple iterators from a vtkDataArray.
typename detail::GetAPITypeImpl< ArrayType, ForceValueTypeForVtkDataArray >::APIType GetAPIType

Note that ArrayType is generic in the above function. When vtk::DataArrayTupleRange is given a vtkDataArray pointer, the generated code produces iterators and reference proxies that rely on the vtkDataArray API. However, when a more derived ArrayType is passed in (for example, vtkFloatArray), specialized implementations are used that generate highly optimized code.

Performance can be further improved when the number of components in the array is known. By passing a compile-time-constant integer as a template parameter, e.g. vtk::DataArrayTupleRange<3>(array), specializations are enabled that allow the compiler to perform additional optimizations.

vtk::DataArrayTupleRange takes an additional two arguments that can be used to restrict the range of tuples to [start, end).

There is a compiler definition / CMake option called VTK_DEBUG_RANGE_ITERATORS that enables checks for proper usage of the range/iterator/reference classes. This slows things down significantly, but is useful for diagnosing problems.

In some situations, developers may want to build in Debug mode while still maintaining decent performance for data-heavy computations. For these usecases, an additional CMake option VTK_ALWAYS_OPTIMIZE_ARRAY_ITERATORS may be enabled to force optimization of code using these iterators. This option will force inlining and enable -O3 (or equivalent) optimization level for iterator code when compiling on platforms that support these features. This option has no effect when VTK_DEBUG_RANGE_ITERATORS is enabled.

Warning
Use caution when using auto to hold values or references obtained from iterators, as they may not behave as expected. This is a deficiency in C++ that affects all proxy iterators (such as those from vector<bool>) that use a reference object instead of an actual C++ reference type. When in doubt, use std::iterator_traits (along with decltype) or the typedefs listed below to determine the proper value/reference type to use. The examples below show how these may be used.

To mitigate this, the following types are defined on the range object:

  • Range::TupleIteratorType: Iterator that visits tuples.
  • Range::ConstTupleIteratorType: Const iterator that visits tuples.
  • Range::TupleReferenceType: Mutable tuple proxy reference.
  • Range::ConstTupleReferenceType: Const tuple proxy reference.
  • Range::ComponentIteratorType: Iterator that visits components in a tuple.
  • Range::ConstComponentIteratorType: Const iterator that visits tuple components.
  • Range::ComponentReferenceType: Reference proxy to a single tuple component.
  • Range::ConstComponentReferenceType: Const reference proxy to a single tuple component.
  • Range::ComponentType: ValueType of components.

These can be accessed via the range objects, e.g.:

auto range = vtk::DataArrayTupleRange(array);
using TupleRef = typename decltype(range)::TupleReferenceType;
using ComponentRef = typename decltype(range)::ComponentReferenceType;
for (TupleRef tuple : range)
{
for (ComponentRef comp : tuple)
{
comp = comp - 1; // Array is modified.
}
}
using ConstTupleRef = typename decltype(range)::ConstTupleReferenceType;
using ComponentType = typename decltype(range)::ComponentType;
for (ConstTupleRef tuple : range)
{
for (ComponentType comp : tuple)
{
comp = comp - 1; // Array is not modified.
}
}
@ range
Definition vtkX3D.h:238
Todo:
Just like the DataArrayValueRange, the tuple range can also accept a forced value type for generic vtkDataArray.

Definition at line 251 of file vtkDataArrayRange.h.

◆ DataArrayValueRange()

template<ComponentIdType TupleSize = detail::DynamicTupleSize, typename ForceValueTypeForVtkDataArray = double, typename ArrayTypePtr = vtkDataArray*>
VTK_ITER_INLINE auto vtk::DataArrayValueRange ( const ArrayTypePtr &  array,
ValueIdType  start = -1,
ValueIdType  end = -1 
) -> typename detail::SelectValueRange<ArrayTypePtr, TupleSize, ForceValueTypeForVtkDataArray>::type

Generate an stl and for-range compatible range of flat AOS iterators from a vtkDataArray.

This function returns a ValueRange object that is compatible with C++11 for-range syntax. The array is traversed as if calling vtkGenericDataArray::GetValue with consecutive, increasing indices. As an example usage, consider a function that takes some instance of vtkDataArray (or a subclass) and sums the values it contains:

template <typename ArrayType>
auto ComputeSum(ArrayType *array) -> vtk::GetAPIType<ArrayType>
{
T sum = 0.;
for (const T val : vtk::DataArrayValueRange(array))
{
sum += val;
}
return sum;
}
VTK_ITER_INLINE auto DataArrayValueRange(const ArrayTypePtr &array, ValueIdType start=-1, ValueIdType end=-1) -> typename detail::SelectValueRange< ArrayTypePtr, TupleSize, ForceValueTypeForVtkDataArray >::type
Generate an stl and for-range compatible range of flat AOS iterators from a vtkDataArray.

These ranges may also be used with STL algorithms:

template <typename ArrayType>
auto ComputeSum(ArrayType *array) -> vtk::GetAPIType<ArrayType>
{
const auto range = vtk::DataArrayValueRange(array);
return std::accumulate(range.begin(), range.end(), 0);
}

Note that ArrayType is generic in the above function. When vtk::DataArrayValueRange is given a vtkDataArray pointer, the generated code produces iterators and reference proxies that rely on the vtkDataArray API. However, when a more derived ArrayType is passed in (for example, vtkFloatArray), specialized implementations are used that generate highly optimized code.

Performance can be further improved when the number of components in the array is known. By passing a compile-time-constant integer as a template parameter, e.g. vtk::DataArrayValueRange<3>(array), specializations are enabled that allow the compiler to perform additional optimizations.

vtk::DataArrayValueRange takes an additional two arguments that can be used to restrict the range of values to [start, end).

There is a compiler definition / CMake option called VTK_DEBUG_RANGE_ITERATORS that enables checks for proper usage of the range/iterator/reference classes. This slows things down significantly, but is useful for diagnosing problems.

In some situations, developers may want to build in Debug mode while still maintaining decent performance for data-heavy computations. For these usecases, an additional CMake option VTK_ALWAYS_OPTIMIZE_ARRAY_ITERATORS may be enabled to force optimization of code using these iterators. This option will force inlining and enable -O3 (or equivalent) optimization level for iterator code when compiling on platforms that support these features. This option has no effect when VTK_DEBUG_RANGE_ITERATORS is enabled.

Warning
Use caution when using auto to hold values or references obtained from iterators, as they may not behave as expected. This is a deficiency in C++ that affects all proxy iterators (such as those from vector<bool>) that use a reference object instead of an actual C++ reference type. When in doubt, use std::iterator_traits (along with decltype) or the typedefs listed below to determine the proper value/reference type to use. The examples below show how these may be used.

To mitigate this, the following types are defined on the range object:

  • Range::IteratorType: Iterator that visits values in AOS order.
  • Range::ConstIteratorType: Const iterator that visits values in AOS order.
  • Range::ReferenceType: Mutable value proxy reference.
  • Range::ConstReferenceType: Const value proxy reference.
  • Range::ValueType: ValueType of array's API.

These can be accessed via the range objects, e.g.:

auto range = vtk::DataArrayValueRange(array);
using RefType = typename decltype(range)::ReferenceType;
for (RefType ref : range)
{ // `ref` is a reference (or reference proxy) to the data held by the array.
ref -= 1; // Array is modified.
}
using ValueType = typename decltype(range)::ValueType;
for (ValueType value : range)
{ // implicitly converts from a reference (or proxy) to a local lvalue `value`
value -= 1; // Array is not modified.
}

Definition at line 361 of file vtkDataArrayRange.h.

◆ Inherits() [1/2]

template<typename VTKObjectType , typename Container >
void vtk::Inherits ( Container &  container)

Populate the container with the name of this class and its ancestors.

The VTKObjectType template-parameter should be a subclass of vtkObjectBase that uses the vtkTypeMacro() to define a Superclass type-alias, as this is how the inheritance hierarchy is traversed.

The version of this function that accepts 2 template parameters uses the second parameter to iterate over a partial hierarchy truncated at (not including) the StopAtType.

Definition at line 170 of file vtkInherits.h.

◆ Inherits() [2/2]

template<typename VTKObjectType , typename StopAtType , typename Container >
void vtk::Inherits ( Container &  container)

Populate the container with the name of this class and its ancestors.

The VTKObjectType template-parameter should be a subclass of vtkObjectBase that uses the vtkTypeMacro() to define a Superclass type-alias, as this is how the inheritance hierarchy is traversed.

The version of this function that accepts 2 template parameters uses the second parameter to iterate over a partial hierarchy truncated at (not including) the StopAtType.

Definition at line 177 of file vtkInherits.h.

◆ Range()

template<typename IterablePtr , typename... Options>
auto vtk::Range ( IterablePtr  iterable,
Options &&...  opts 
) -> typename detail::IterableTraits<typename detail::StripPointers<IterablePtr>::type>::RangeType

Generate an iterable STL proxy object for a VTK container.

Currently supports:

Usage:

for (auto item : vtk::Range(myCollection))
{
// Use item.
}
// or:
using Opts = vtk::vtkDataObjectTreeOptions;
auto range = vtk::Range(dataObjTree,
Opts::TraverseSubTree | Opts::VisitOnlyLeaves);
some_algo(range.begin(), range.end());
auto Range(IterablePtr iterable, Options &&... opts) -> typename detail::IterableTraits< typename detail::StripPointers< IterablePtr >::type >::RangeType
Generate an iterable STL proxy object for a VTK container.
Definition vtkRange.h:74

Definition at line 74 of file vtkRange.h.

◆ MakeSmartPointer()

template<typename T >
vtkSmartPointer< T > vtk::MakeSmartPointer ( T *  obj)

Construct a vtkSmartPointer<T> containing obj.

A new reference is added to obj.

Definition at line 468 of file vtkSmartPointer.h.

◆ TakeSmartPointer()

template<typename T >
vtkSmartPointer< T > vtk::TakeSmartPointer ( T *  obj)

Construct a vtkSmartPointer<T> containing obj.

obj's reference count is not changed.

Definition at line 476 of file vtkSmartPointer.h.

◆ CreateStructuredPointArray()

template<typename ValueType >
vtkSmartPointer< vtkImplicitArray< vtkStructuredPointBackend< ValueType > > > vtk::CreateStructuredPointArray ( vtkDataArray xCoords,
vtkDataArray yCoords,
vtkDataArray zCoords,
int  extent[6],
int  dataDescription,
double  dirMatrix[9] 
)

Create an implicit point array from the given coordinates and direction matrix which is optional.

xCoords, yCoords and zCoords are the coordinates of the points. extent is the extent of the dataset. dataDescription is the data description of the dataset. dirMatrix is the direction matrix of the dataset (if any, else provide a homogeneous matrix).

◆ TypeName()

template<typename ObjectType >
std::string vtk::TypeName ( )
inline

Return the demangled type-name of the provided ObjectType.

Note that if the <cxxabi.h> header is not present or does not provide abi::__cxa_demangle(), a mangled name will be returned.

Definition at line 116 of file vtkTypeName.h.

◆ TypeToken()

template<typename ObjectType >
vtkStringToken vtk::TypeToken ( )
inline

Return a string token holding a hash of the demangled type-name of the provided ObjectType.

Note that this function should become constexpr, so the hash will be computed at compile time (once VTK allows c++14 extensions). Because of this, the string for the hash cannot be added to the vtkStringManager and thus may cause exceptions later.

Definition at line 132 of file vtkTypeName.h.

◆ TakeWeakPointer()

template<typename T >
vtkWeakPointer< T > vtk::TakeWeakPointer ( T *  obj)

Construct a vtkWeakPointer<T> containing obj.

obj's reference count is not changed.

Definition at line 295 of file vtkWeakPointer.h.