VTK  9.1.0
DataArrayConverters.h
Go to the documentation of this file.
1//=============================================================================
2//
3// Copyright (c) Kitware, Inc.
4// All rights reserved.
5// See LICENSE.txt for details.
6//
7// This software is distributed WITHOUT ANY WARRANTY; without even
8// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
9// PURPOSE. See the above copyright notice for more information.
10//
11// Copyright 2012 Sandia Corporation.
12// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
13// the U.S. Government retains certain rights in this software.
14//
15//=============================================================================
16
17#ifndef vtkmlib_DataArrayConverters_h
18#define vtkmlib_DataArrayConverters_h
19
20#include "vtkAcceleratorsVTKmCoreModule.h" //required for correct implementation
21#include "vtkmConfigCore.h" //required for general vtkm setup
22
25
26#include <vtkm/cont/ArrayHandleSOA.h>
27#include <vtkm/cont/Field.h>
28#include <vtkm/cont/UnknownArrayHandle.h>
29
30#include <type_traits> // for std::underlying_type
31
32class vtkDataArray;
33class vtkPoints;
34
35namespace vtkm
36{
37namespace cont
38{
39class CoordinateSystem;
40}
41}
42
43namespace tovtkm
44{
45
46template <typename DataArrayType, vtkm::IdComponent NumComponents>
48
49template <typename T, vtkm::IdComponent NumComponents>
51{
52 using ValueType =
53 typename std::conditional<NumComponents == 1, T, vtkm::Vec<T, NumComponents>>::type;
54 using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagBasic>;
55 using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagBasic>;
56
58 {
59 return vtkm::cont::make_ArrayHandle(reinterpret_cast<ValueType*>(input->GetPointer(0)),
60 input->GetNumberOfTuples(), vtkm::CopyFlag::Off);
61 }
62};
63
64template <typename T, vtkm::IdComponent NumComponents>
66{
67 using ValueType = vtkm::Vec<T, NumComponents>;
68 using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagSOA>;
69 using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagSOA>;
70
72 {
73 vtkm::Id numValues = input->GetNumberOfTuples();
74 vtkm::cont::ArrayHandleSOA<ValueType> handle;
75 for (vtkm::IdComponent i = 0; i < NumComponents; ++i)
76 {
77 handle.SetArray(i,
78 vtkm::cont::make_ArrayHandle<T>(reinterpret_cast<T*>(input->GetComponentArrayPointer(i)),
79 numValues, vtkm::CopyFlag::Off));
80 }
81
82 return handle;
83 }
84};
85
86template <typename T>
88{
89 using StorageType = vtkm::cont::internal::Storage<T, vtkm::cont::StorageTagBasic>;
90 using ArrayHandleType = vtkm::cont::ArrayHandle<T, vtkm::cont::StorageTagBasic>;
91
93 {
94 return vtkm::cont::make_ArrayHandle(
95 input->GetComponentArrayPointer(0), input->GetNumberOfTuples(), vtkm::CopyFlag::Off);
96 }
97};
98
99enum class FieldsFlag
100{
101 None = 0x0,
102 Points = 0x1,
103 Cells = 0x2,
104
106};
107
108}
109
110namespace fromvtkm
111{
112
113VTKACCELERATORSVTKMCORE_EXPORT
114vtkDataArray* Convert(const vtkm::cont::Field& input);
115
116VTKACCELERATORSVTKMCORE_EXPORT
117vtkDataArray* Convert(const vtkm::cont::UnknownArrayHandle& input, const char* name);
118
119VTKACCELERATORSVTKMCORE_EXPORT
120vtkPoints* Convert(const vtkm::cont::CoordinateSystem& input);
121
122}
123
125{
126 using T = std::underlying_type<tovtkm::FieldsFlag>::type;
127 return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) & static_cast<T>(b));
128}
129
131{
132 using T = std::underlying_type<tovtkm::FieldsFlag>::type;
133 return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) | static_cast<T>(b));
134}
135
136#endif // vtkmlib_ArrayConverters_h
tovtkm::FieldsFlag operator&(const tovtkm::FieldsFlag &a, const tovtkm::FieldsFlag &b)
tovtkm::FieldsFlag operator|(const tovtkm::FieldsFlag &a, const tovtkm::FieldsFlag &b)
Array-Of-Structs implementation of vtkGenericDataArray.
ValueType * GetPointer(vtkIdType valueIdx)
Get the address of a particular data index.
vtkIdType GetNumberOfTuples() const
Get the number of complete tuples (a component group) in the array.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
represent and manipulate 3D points
Definition: vtkPoints.h:143
Struct-Of-Arrays implementation of vtkGenericDataArray.
ValueType * GetComponentArrayPointer(int comp)
Return a pointer to a contiguous block of memory containing all values for a particular components (i...
VTKACCELERATORSVTKMCORE_EXPORT vtkDataArray * Convert(const vtkm::cont::Field &input)
typename std::conditional< NumComponents==1, T, vtkm::Vec< T, NumComponents > >::type ValueType
vtkm::cont::ArrayHandle< ValueType, vtkm::cont::StorageTagBasic > ArrayHandleType
vtkm::cont::internal::Storage< ValueType, vtkm::cont::StorageTagBasic > StorageType
static ArrayHandleType Wrap(vtkAOSDataArrayTemplate< T > *input)
static ArrayHandleType Wrap(vtkSOADataArrayTemplate< T > *input)
vtkm::cont::ArrayHandle< T, vtkm::cont::StorageTagBasic > ArrayHandleType
vtkm::cont::internal::Storage< T, vtkm::cont::StorageTagBasic > StorageType
vtkm::cont::internal::Storage< ValueType, vtkm::cont::StorageTagSOA > StorageType
static ArrayHandleType Wrap(vtkSOADataArrayTemplate< T > *input)
vtkm::cont::ArrayHandle< ValueType, vtkm::cont::StorageTagSOA > ArrayHandleType