VTK  9.2.6
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
49inline static const char* NoNameVTKFieldName()
50{
51 static const char* name = "NoNameVTKField";
52 return name;
53}
54
55template <typename DataArrayType, vtkm::IdComponent NumComponents>
57
58template <typename T, vtkm::IdComponent NumComponents>
60{
61 using ValueType =
62 typename std::conditional<NumComponents == 1, T, vtkm::Vec<T, NumComponents>>::type;
63 using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagBasic>;
64 using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagBasic>;
65
67 {
68 return vtkm::cont::make_ArrayHandle(reinterpret_cast<ValueType*>(input->GetPointer(0)),
69 input->GetNumberOfTuples(), vtkm::CopyFlag::Off);
70 }
71};
72
73template <typename T, vtkm::IdComponent NumComponents>
75{
76 using ValueType = vtkm::Vec<T, NumComponents>;
77 using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagSOA>;
78 using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagSOA>;
79
81 {
82 vtkm::Id numValues = input->GetNumberOfTuples();
83 vtkm::cont::ArrayHandleSOA<ValueType> handle;
84 for (vtkm::IdComponent i = 0; i < NumComponents; ++i)
85 {
86 handle.SetArray(i,
87 vtkm::cont::make_ArrayHandle<T>(reinterpret_cast<T*>(input->GetComponentArrayPointer(i)),
88 numValues, vtkm::CopyFlag::Off));
89 }
90
91 return handle;
92 }
93};
94
95template <typename T>
97{
98 using StorageType = vtkm::cont::internal::Storage<T, vtkm::cont::StorageTagBasic>;
99 using ArrayHandleType = vtkm::cont::ArrayHandle<T, vtkm::cont::StorageTagBasic>;
100
102 {
103 return vtkm::cont::make_ArrayHandle(
104 input->GetComponentArrayPointer(0), input->GetNumberOfTuples(), vtkm::CopyFlag::Off);
105 }
106};
107
108enum class FieldsFlag
109{
110 None = 0x0,
111 Points = 0x1,
112 Cells = 0x2,
113
115};
116
117}
118
119namespace fromvtkm
120{
121
122VTKACCELERATORSVTKMCORE_EXPORT
123vtkDataArray* Convert(const vtkm::cont::Field& input);
124
125VTKACCELERATORSVTKMCORE_EXPORT
126vtkDataArray* Convert(const vtkm::cont::UnknownArrayHandle& input, const char* name);
127
128VTKACCELERATORSVTKMCORE_EXPORT
129vtkPoints* Convert(const vtkm::cont::CoordinateSystem& input);
130
131}
132
134{
135 using T = std::underlying_type<tovtkm::FieldsFlag>::type;
136 return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) & static_cast<T>(b));
137}
138
140{
141 using T = std::underlying_type<tovtkm::FieldsFlag>::type;
142 return static_cast<tovtkm::FieldsFlag>(static_cast<T>(a) | static_cast<T>(b));
143}
144
145#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
represent and manipulate 3D points
Definition vtkPoints.h:40
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)
static const char * NoNameVTKFieldName()
Temporary name for arrays converted from VTK that do not have a name.
vtkm::cont::ArrayHandle< ValueType, vtkm::cont::StorageTagBasic > ArrayHandleType
typename std::conditional< NumComponents==1, T, vtkm::Vec< T, NumComponents > >::type ValueType
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