VTK  9.3.1
DataArrayConverters.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) Kitware, Inc.
3 // SPDX-FileCopyrightText: Copyright 2012 Sandia Corporation.
4 // SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-USGov
5 
6 #ifndef vtkmlib_DataArrayConverters_h
7 #define vtkmlib_DataArrayConverters_h
8 
9 #include "vtkAcceleratorsVTKmCoreModule.h" //required for correct implementation
10 #include "vtkmConfigCore.h" //required for general vtkm setup
11 
14 
15 #include <vtkm/cont/ArrayHandleSOA.h>
16 #include <vtkm/cont/Field.h>
17 #include <vtkm/cont/UnknownArrayHandle.h>
18 
19 #include <type_traits> // for std::underlying_type
20 
21 namespace vtkm
22 {
23 namespace cont
24 {
25 class CoordinateSystem;
26 }
27 }
28 
29 VTK_ABI_NAMESPACE_BEGIN
30 class vtkDataArray;
31 class vtkPoints;
32 VTK_ABI_NAMESPACE_END
33 
34 namespace tovtkm
35 {
36 VTK_ABI_NAMESPACE_BEGIN
37 
41 inline static const char* NoNameVTKFieldName()
42 {
43  static const char* name = "NoNameVTKField";
44  return name;
45 }
46 
47 template <typename DataArrayType, vtkm::IdComponent NumComponents>
49 
50 template <typename T, vtkm::IdComponent NumComponents>
52 {
53  using ValueType =
54  typename std::conditional<NumComponents == 1, T, vtkm::Vec<T, NumComponents>>::type;
55  using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagBasic>;
56  using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagBasic>;
57 
59  {
60  return vtkm::cont::make_ArrayHandle(reinterpret_cast<ValueType*>(input->GetPointer(0)),
61  input->GetNumberOfTuples(), vtkm::CopyFlag::Off);
62  }
63 };
64 
65 template <typename T, vtkm::IdComponent NumComponents>
67 {
68  using ValueType = vtkm::Vec<T, NumComponents>;
69  using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagSOA>;
70  using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagSOA>;
71 
73  {
74  vtkm::Id numValues = input->GetNumberOfTuples();
75  vtkm::cont::ArrayHandleSOA<ValueType> handle;
76  for (vtkm::IdComponent i = 0; i < NumComponents; ++i)
77  {
78  handle.SetArray(i,
79  vtkm::cont::make_ArrayHandle<T>(reinterpret_cast<T*>(input->GetComponentArrayPointer(i)),
80  numValues, vtkm::CopyFlag::Off));
81  }
82 
83  return std::move(handle);
84  }
85 };
86 
87 template <typename T>
89 {
90  using StorageType = vtkm::cont::internal::Storage<T, vtkm::cont::StorageTagBasic>;
91  using ArrayHandleType = vtkm::cont::ArrayHandle<T, vtkm::cont::StorageTagBasic>;
92 
94  {
95  return vtkm::cont::make_ArrayHandle(
96  input->GetComponentArrayPointer(0), input->GetNumberOfTuples(), vtkm::CopyFlag::Off);
97  }
98 };
99 
100 enum class FieldsFlag
101 {
102  None = 0x0,
103  Points = 0x1,
104  Cells = 0x2,
105 
107 };
108 
109 VTK_ABI_NAMESPACE_END
110 }
111 
112 namespace fromvtkm
113 {
114 VTK_ABI_NAMESPACE_BEGIN
115 
116 VTKACCELERATORSVTKMCORE_EXPORT
117 vtkDataArray* Convert(const vtkm::cont::Field& input);
118 
119 VTKACCELERATORSVTKMCORE_EXPORT
120 vtkDataArray* Convert(const vtkm::cont::UnknownArrayHandle& input, const char* name);
121 
122 VTKACCELERATORSVTKMCORE_EXPORT
123 vtkPoints* Convert(const vtkm::cont::CoordinateSystem& input);
124 
125 VTK_ABI_NAMESPACE_END
126 }
127 
128 VTK_ABI_NAMESPACE_BEGIN
130 {
132  return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) & static_cast<T>(b));
133 }
134 
136 {
138  return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) | static_cast<T>(b));
139 }
140 VTK_ABI_NAMESPACE_END
141 
142 #endif // vtkmlib_ArrayConverters_h
143 /* VTK-HeaderTest-Exclude: DataArrayConverters.h */
Struct-Of-Arrays implementation of vtkGenericDataArray.
static ArrayHandleType Wrap(vtkSOADataArrayTemplate< T > *input)
tovtkm::FieldsFlag operator&(const tovtkm::FieldsFlag &a, const tovtkm::FieldsFlag &b)
ValueType * GetPointer(vtkIdType valueIdx)
Get the address of a particular data index.
vtkm::cont::internal::Storage< ValueType, vtkm::cont::StorageTagBasic > StorageType
tovtkm::FieldsFlag operator|(const tovtkm::FieldsFlag &a, const tovtkm::FieldsFlag &b)
static ArrayHandleType Wrap(vtkSOADataArrayTemplate< T > *input)
vtkm::cont::internal::Storage< ValueType, vtkm::cont::StorageTagSOA > StorageType
ValueType * GetComponentArrayPointer(int comp)
Return a pointer to a contiguous block of memory containing all values for a particular components (i...
static const char * NoNameVTKFieldName()
Temporary name for arrays converted from VTK that do not have a name.
vtkm::cont::ArrayHandle< T, vtkm::cont::StorageTagBasic > ArrayHandleType
vtkm::cont::ArrayHandle< ValueType, vtkm::cont::StorageTagSOA > ArrayHandleType
Array-Of-Structs implementation of vtkGenericDataArray.
VTKACCELERATORSVTKMCORE_EXPORT vtkDataArray * Convert(const vtkm::cont::Field &input)
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
vtkm::cont::ArrayHandle< ValueType, vtkm::cont::StorageTagBasic > ArrayHandleType
vtkm::cont::internal::Storage< T, vtkm::cont::StorageTagBasic > StorageType
vtkIdType GetNumberOfTuples() const
Get the number of complete tuples (a component group) in the array.
static ArrayHandleType Wrap(vtkAOSDataArrayTemplate< T > *input)
typename std::conditional< NumComponents==1, T, vtkm::Vec< T, NumComponents >>::type ValueType
represent and manipulate 3D points
Definition: vtkPoints.h:28