VTK  9.0.1
ArrayConverters.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_ArrayConverters_h
18 #define vtkmlib_ArrayConverters_h
19 
20 #include "vtkAcceleratorsVTKmModule.h" //required for correct implementation
21 #include "vtkmConfig.h" //required for general vtkm setup
22 
25 
26 #include <vtkm/cont/ArrayHandleSOA.h>
27 #include <vtkm/cont/Field.h>
28 
29 #include <type_traits> // for std::underlying_type
30 
31 class vtkDataArray;
32 class vtkDataSet;
33 class vtkPoints;
34 
35 namespace vtkm
36 {
37 namespace cont
38 {
39 class DataSet;
40 class CoordinateSystem;
41 }
42 }
43 
44 namespace tovtkm
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(
61  reinterpret_cast<ValueType*>(input->GetPointer(0)), input->GetNumberOfTuples());
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::internal::Storage<ValueType, vtkm::cont::StorageTagSOA> storage;
76  for (vtkm::IdComponent i = 0; i < NumComponents; ++i)
77  {
78  storage.SetArray(i,
79  vtkm::cont::make_ArrayHandle<T>(
80  reinterpret_cast<T*>(input->GetComponentArrayPointer(i)), numValues));
81  }
82 
83  return vtkm::cont::ArrayHandleSOA<ValueType>(storage);
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());
97  }
98 };
99 
100 enum class FieldsFlag
101 {
102  None = 0x0,
103  Points = 0x1,
104  Cells = 0x2,
105 
107 };
108 
109 VTKACCELERATORSVTKM_EXPORT
111 
112 // determine the type and call the proper Convert routine
113 VTKACCELERATORSVTKM_EXPORT
114 vtkm::cont::Field Convert(vtkDataArray* input, int association);
115 }
116 
117 namespace fromvtkm
118 {
119 
120 VTKACCELERATORSVTKM_EXPORT
121 vtkDataArray* Convert(const vtkm::cont::Field& input);
122 
123 VTKACCELERATORSVTKM_EXPORT
124 vtkPoints* Convert(const vtkm::cont::CoordinateSystem& input);
125 
126 VTKACCELERATORSVTKM_EXPORT
127 bool ConvertArrays(const vtkm::cont::DataSet& input, vtkDataSet* output);
128 }
129 
131 {
133  return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) & static_cast<T>(b));
134 }
135 
137 {
139  return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) | static_cast<T>(b));
140 }
141 
142 #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:50
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
represent and manipulate 3D points
Definition: vtkPoints.h:34
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...
VTKACCELERATORSVTKM_EXPORT bool ConvertArrays(const vtkm::cont::DataSet &input, vtkDataSet *output)
VTKACCELERATORSVTKM_EXPORT vtkDataArray * Convert(const vtkm::cont::Field &input)
VTKACCELERATORSVTKM_EXPORT vtkm::cont::Field Convert(vtkDataArray *input, int association)
VTKACCELERATORSVTKM_EXPORT void ProcessFields(vtkDataSet *input, vtkm::cont::DataSet &dataset, tovtkm::FieldsFlag fields)
@ type
Definition: vtkX3D.h:522
std::map< std::string, DataArray > DataSet
key: variable name, value: DataArray
Definition: VTXTypes.h:39
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