VTK
9.3.1
|
A backend for the vtkImplicitArray
framework allowing one to use a subset of a given data array, by providing a vtkIdList
or vtkDataArray
of indexes as indirection, as another vtkDataArray
without any excess memory consumption.
More...
#include <vtkImplicitArray.h>
Public Member Functions | |
~vtkIndexedImplicitBackend () | |
ValueType | operator() (int idx) const |
Indexing operation for the indexed array respecting the backend expectations of vtkImplicitArray More... | |
vtkIndexedImplicitBackend (vtkIdList *indexes, vtkDataArray *array) | |
Constructor. More... | |
vtkIndexedImplicitBackend (vtkDataArray *indexes, vtkDataArray *array) | |
Constructor. More... | |
A backend for the vtkImplicitArray
framework allowing one to use a subset of a given data array, by providing a vtkIdList
or vtkDataArray
of indexes as indirection, as another vtkDataArray
without any excess memory consumption.
This structure can be classified as a closure and can be called using syntax similar to a function call.
An example of potential usage in a vtkImplicitArray
:
``` vtkNew<vtkIntArray> baseArray; baseArray->SetNumberOfComponents(1); baseArray->SetNumberOfTuples(100); auto range = vtk::DataArrayValueRange<1>(baseArray); std::iota(range.begin(), range.end(), 0);
vtkNew<vtkIdList> handles; handles->SetNumberOfIds(100); for (vtkIdType idx = 0; idx < 100; idx++) { handles->SetId(idx, 99-idx); }
vtkNew<vtkImplicitArray<vtkIndexedImplicitBackend<int>>> indexedArr; // More compact with // vtkIndexedArray<int>
// if available indexedArr->SetBackend(std::make_shared<vtkIndexedImplicitBackend<int>>(handles, baseArray)); indexedArr->SetNumberOfComponents(1); indexedArr->SetNumberOfTuples(100); CHECK(indexedArr->GetValue(57) == 42);
```
Definition at line 508 of file vtkImplicitArray.h.
vtkIndexedImplicitBackend< ValueType >::vtkIndexedImplicitBackend | ( | vtkIdList * | indexes, |
vtkDataArray * | array | ||
) |
Constructor.
indexes | list of indexes to use for indirection of the array |
array | base array of interest |
vtkIndexedImplicitBackend< ValueType >::vtkIndexedImplicitBackend | ( | vtkDataArray * | indexes, |
vtkDataArray * | array | ||
) |
Constructor.
indexes | list of indexes to use for indirection of the array |
array | base array of interest |
vtkIndexedImplicitBackend< ValueType >::~vtkIndexedImplicitBackend | ( | ) |
ValueType vtkIndexedImplicitBackend< ValueType >::operator() | ( | int | idx | ) | const |
Indexing operation for the indexed array respecting the backend expectations of vtkImplicitArray