VTK-m  2.2
DataSetBuilderExplicit.h
Go to the documentation of this file.
1 //============================================================================
2 // Copyright (c) Kitware, Inc.
3 // All rights reserved.
4 // See LICENSE.txt for details.
5 //
6 // This software is distributed WITHOUT ANY WARRANTY; without even
7 // the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
8 // PURPOSE. See the above copyright notice for more information.
9 //============================================================================
10 #ifndef vtk_m_cont_DataSetBuilderExplicit_h
11 #define vtk_m_cont_DataSetBuilderExplicit_h
12 
18 #include <vtkm/cont/DataSet.h>
19 
20 namespace vtkm
21 {
22 namespace cont
23 {
24 
25 //Coordinates builder??
26 
28 {
29 public:
30  VTKM_CONT
32 
51  template <typename T>
52  VTKM_CONT static vtkm::cont::DataSet Create(const std::vector<T>& xVals,
53  const std::vector<vtkm::UInt8>& shapes,
54  const std::vector<vtkm::IdComponent>& numIndices,
55  const std::vector<vtkm::Id>& connectivity,
56  const std::string& coordsNm = "coords")
57  {
58  std::vector<T> yVals(xVals.size(), 0), zVals(xVals.size(), 0);
60  xVals, yVals, zVals, shapes, numIndices, connectivity, coordsNm);
61  }
62 
82  template <typename T>
83  VTKM_CONT static vtkm::cont::DataSet Create(const std::vector<T>& xVals,
84  const std::vector<T>& yVals,
85  const std::vector<vtkm::UInt8>& shapes,
86  const std::vector<vtkm::IdComponent>& numIndices,
87  const std::vector<vtkm::Id>& connectivity,
88  const std::string& coordsNm = "coords")
89  {
90  std::vector<T> zVals(xVals.size(), 0);
92  xVals, yVals, zVals, shapes, numIndices, connectivity, coordsNm);
93  }
94 
115  template <typename T>
116  VTKM_CONT static vtkm::cont::DataSet Create(const std::vector<T>& xVals,
117  const std::vector<T>& yVals,
118  const std::vector<T>& zVals,
119  const std::vector<vtkm::UInt8>& shapes,
120  const std::vector<vtkm::IdComponent>& numIndices,
121  const std::vector<vtkm::Id>& connectivity,
122  const std::string& coordsNm = "coords");
123 
142  template <typename T>
143  VTKM_CONT static vtkm::cont::DataSet Create(const std::vector<vtkm::Vec<T, 3>>& coords,
144  const std::vector<vtkm::UInt8>& shapes,
145  const std::vector<vtkm::IdComponent>& numIndices,
146  const std::vector<vtkm::Id>& connectivity,
147  const std::string& coordsNm = "coords");
148 
168  template <typename T>
173  const vtkm::cont::ArrayHandle<vtkm::Id>& connectivity,
174  const std::string& coordsNm = "coords")
175  {
176  auto offsets = vtkm::cont::ConvertNumComponentsToOffsets(numIndices);
177  return DataSetBuilderExplicit::BuildDataSet(coords, shapes, offsets, connectivity, coordsNm);
178  }
179 
201  template <typename T, typename CellShapeTag>
202  VTKM_CONT static vtkm::cont::DataSet Create(const std::vector<vtkm::Vec<T, 3>>& coords,
203  CellShapeTag tag,
204  vtkm::IdComponent numberOfPointsPerCell,
205  const std::vector<vtkm::Id>& connectivity,
206  const std::string& coordsNm = "coords");
207 
229  template <typename T, typename CellShapeTag>
232  CellShapeTag tag,
233  vtkm::IdComponent numberOfPointsPerCell,
234  const vtkm::cont::ArrayHandle<vtkm::Id>& connectivity,
235  const std::string& coordsNm = "coords")
236  {
238  coords, tag, numberOfPointsPerCell, connectivity, coordsNm);
239  }
240 
241 private:
242  template <typename T>
243  VTKM_CONT static vtkm::cont::DataSet BuildDataSet(
246  const vtkm::cont::ArrayHandle<vtkm::Id>& offsets,
247  const vtkm::cont::ArrayHandle<vtkm::Id>& connectivity,
248  const std::string& coordsNm);
249 
250  template <typename T, typename CellShapeTag>
251  VTKM_CONT static vtkm::cont::DataSet BuildDataSet(
253  CellShapeTag tag,
254  vtkm::IdComponent numberOfPointsPerCell,
255  const vtkm::cont::ArrayHandle<vtkm::Id>& connectivity,
256  const std::string& coordsNm);
257 };
258 
259 template <typename T>
261  const std::vector<T>& xVals,
262  const std::vector<T>& yVals,
263  const std::vector<T>& zVals,
264  const std::vector<vtkm::UInt8>& shapes,
265  const std::vector<vtkm::IdComponent>& numIndices,
266  const std::vector<vtkm::Id>& connectivity,
267  const std::string& coordsNm)
268 {
269  VTKM_ASSERT(xVals.size() == yVals.size() && yVals.size() == zVals.size() && xVals.size() > 0);
270 
272  coordsArray.Allocate(static_cast<vtkm::Id>(xVals.size()));
273  auto coordsPortal = coordsArray.WritePortal();
274  for (std::size_t index = 0; index < xVals.size(); ++index)
275  {
276  coordsPortal.Set(static_cast<vtkm::Id>(index),
277  vtkm::make_Vec(static_cast<vtkm::FloatDefault>(xVals[index]),
278  static_cast<vtkm::FloatDefault>(yVals[index]),
279  static_cast<vtkm::FloatDefault>(zVals[index])));
280  }
281 
282  auto shapesArray = vtkm::cont::make_ArrayHandle(shapes, vtkm::CopyFlag::On);
283  auto connArray = vtkm::cont::make_ArrayHandle(connectivity, vtkm::CopyFlag::On);
284 
285  auto offsetsArray = vtkm::cont::ConvertNumComponentsToOffsets(
287 
289  coordsArray, shapesArray, offsetsArray, connArray, coordsNm);
290 }
291 
292 template <typename T>
294  const std::vector<vtkm::Vec<T, 3>>& coords,
295  const std::vector<vtkm::UInt8>& shapes,
296  const std::vector<vtkm::IdComponent>& numIndices,
297  const std::vector<vtkm::Id>& connectivity,
298  const std::string& coordsNm)
299 {
300  auto coordsArray = vtkm::cont::make_ArrayHandle(coords, vtkm::CopyFlag::On);
301  auto shapesArray = vtkm::cont::make_ArrayHandle(shapes, vtkm::CopyFlag::On);
302  auto connArray = vtkm::cont::make_ArrayHandle(connectivity, vtkm::CopyFlag::On);
303 
304  auto offsetsArray = vtkm::cont::ConvertNumComponentsToOffsets(
306 
308  coordsArray, shapesArray, offsetsArray, connArray, coordsNm);
309 }
310 
311 template <typename T>
315  const vtkm::cont::ArrayHandle<vtkm::Id>& offsets,
316  const vtkm::cont::ArrayHandle<vtkm::Id>& connectivity,
317  const std::string& coordsNm)
318 {
319  vtkm::cont::DataSet dataSet;
320 
321  dataSet.AddCoordinateSystem(vtkm::cont::CoordinateSystem(coordsNm, coords));
322  vtkm::Id nPts = static_cast<vtkm::Id>(coords.GetNumberOfValues());
324 
325  cellSet.Fill(nPts, shapes, connectivity, offsets);
326  dataSet.SetCellSet(cellSet);
327 
328  return dataSet;
329 }
330 
331 template <typename T, typename CellShapeTag>
333  const std::vector<vtkm::Vec<T, 3>>& coords,
334  CellShapeTag tag,
335  vtkm::IdComponent numberOfPointsPerCell,
336  const std::vector<vtkm::Id>& connectivity,
337  const std::string& coordsNm)
338 {
339  auto coordsArray = vtkm::cont::make_ArrayHandle(coords, vtkm::CopyFlag::On);
340  auto connArray = vtkm::cont::make_ArrayHandle(connectivity, vtkm::CopyFlag::On);
341 
343  coordsArray, tag, numberOfPointsPerCell, connArray, coordsNm);
344 }
345 
346 template <typename T, typename CellShapeTag>
349  CellShapeTag tag,
350  vtkm::IdComponent numberOfPointsPerCell,
351  const vtkm::cont::ArrayHandle<vtkm::Id>& connectivity,
352  const std::string& coordsNm)
353 {
354  (void)tag; //C4100 false positive workaround
355  vtkm::cont::DataSet dataSet;
356 
357  dataSet.AddCoordinateSystem(vtkm::cont::CoordinateSystem(coordsNm, coords));
359 
360  cellSet.Fill(coords.GetNumberOfValues(), tag.Id, numberOfPointsPerCell, connectivity);
361  dataSet.SetCellSet(cellSet);
362 
363  return dataSet;
364 }
365 
370 {
371 public:
373 
380  VTKM_CONT void Begin(const std::string& coordName = "coords");
381 
386  VTKM_CONT vtkm::Id AddPoint(const vtkm::Vec3f& pt);
387 
392  template <typename T>
394  {
395  return AddPoint(static_cast<vtkm::Vec3f>(pt));
396  }
397 
404  VTKM_CONT vtkm::Id AddPoint(const vtkm::FloatDefault& x,
405  const vtkm::FloatDefault& y,
406  const vtkm::FloatDefault& z = 0);
407 
414  template <typename T>
415  VTKM_CONT vtkm::Id AddPoint(const T& x, const T& y, const T& z = 0)
416  {
417  return AddPoint(static_cast<vtkm::FloatDefault>(x),
418  static_cast<vtkm::FloatDefault>(y),
419  static_cast<vtkm::FloatDefault>(z));
420  }
421 
427  VTKM_CONT void AddCell(const vtkm::UInt8& shape, const std::vector<vtkm::Id>& conn);
428 
435  VTKM_CONT void AddCell(const vtkm::UInt8& shape,
436  const vtkm::Id* conn,
437  const vtkm::IdComponent& n);
438 
445  VTKM_CONT void AddCell(vtkm::UInt8 shape);
446 
450  VTKM_CONT void AddCellPoint(vtkm::Id pointIndex);
451 
457 
458 private:
459  std::string coordNm;
460 
461  std::vector<vtkm::Vec3f> points;
462  std::vector<vtkm::UInt8> shapes;
463  std::vector<vtkm::IdComponent> numIdx;
464  std::vector<vtkm::Id> connectivity;
465 };
466 }
467 }
468 
469 #endif //vtk_m_cont_DataSetBuilderExplicit_h
vtkm::cont::ArrayHandle
Manages an array-worth of data.
Definition: ArrayHandle.h:300
vtkm::cont::CellSetSingleType::Fill
void Fill(vtkm::Id numPoints, vtkm::UInt8 shapeId, vtkm::IdComponent numberOfPointsPerCell, const vtkm::cont::ArrayHandle< vtkm::Id, ConnectivityStorageTag > &connectivity)
Set all the cells of the mesh.
Definition: CellSetSingleType.h:200
vtkm::cont::DataSetBuilderExplicit::Create
static vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle< vtkm::Vec< T, 3 >> &coords, CellShapeTag tag, vtkm::IdComponent numberOfPointsPerCell, const vtkm::cont::ArrayHandle< vtkm::Id > &connectivity, const std::string &coordsNm="coords")
Create a 3D DataSet with arbitrary cell connectivity for a single cell type.
Definition: DataSetBuilderExplicit.h:230
vtkm
Groups connected points that have the same field value.
Definition: Atomic.h:19
ArrayHandleCast.h
VTKM_ASSERT
#define VTKM_ASSERT(condition)
Definition: Assert.h:43
vtkm::make_Vec
constexpr vtkm::Vec< T, vtkm::IdComponent(sizeof...(Ts)+1)> make_Vec(T value0, Ts &&... args)
Initializes and returns a Vec containing all the arguments.
Definition: Types.h:1253
vtkm::IdComponent
vtkm::Int32 IdComponent
Base type to use to index small lists.
Definition: Types.h:194
ArrayPortalToIterators.h
vtkm::cont::DataSetBuilderExplicit::Create
static vtkm::cont::DataSet Create(const vtkm::cont::ArrayHandle< vtkm::Vec< T, 3 >> &coords, const vtkm::cont::ArrayHandle< vtkm::UInt8 > &shapes, const vtkm::cont::ArrayHandle< vtkm::IdComponent > &numIndices, const vtkm::cont::ArrayHandle< vtkm::Id > &connectivity, const std::string &coordsNm="coords")
Create a 3D DataSet with arbitrary cell connectivity.
Definition: DataSetBuilderExplicit.h:169
ArrayHandleCompositeVector.h
vtkm::cont::DataSet::SetCellSet
void SetCellSet(const CellSetType &cellSet)
Definition: DataSet.h:396
vtkm::cont::CellSetSingleType
An explicit cell set with all cells of the same shape.
Definition: CastAndCall.h:34
vtkm::cont::DataSet
Contains and manages the geometric data structures that VTK-m operates on.
Definition: DataSet.h:57
CoordinateSystem.h
vtkm::cont::CellSetExplicit::Fill
void Fill(vtkm::Id numPoints, const vtkm::cont::ArrayHandle< vtkm::UInt8, ShapesStorageTag > &cellTypes, const vtkm::cont::ArrayHandle< vtkm::Id, ConnectivityStorageTag > &connectivity, const vtkm::cont::ArrayHandle< vtkm::Id, OffsetsStorageTag > &offsets)
Set all the cells of the mesh.
vtkm::cont::DataSetBuilderExplicitIterative::AddPoint
vtkm::Id AddPoint(const T &x, const T &y, const T &z=0)
Add a point to the DataSet.
Definition: DataSetBuilderExplicit.h:415
vtkm::cont::CoordinateSystem
Manages a coordinate system for a DataSet.
Definition: CoordinateSystem.h:30
vtkm::cont::DataSetBuilderExplicitIterative::points
std::vector< vtkm::Vec3f > points
Definition: DataSetBuilderExplicit.h:461
vtkm::cont::DataSetBuilderExplicit::BuildDataSet
static vtkm::cont::DataSet BuildDataSet(const vtkm::cont::ArrayHandle< vtkm::Vec< T, 3 >> &coords, const vtkm::cont::ArrayHandle< vtkm::UInt8 > &shapes, const vtkm::cont::ArrayHandle< vtkm::Id > &offsets, const vtkm::cont::ArrayHandle< vtkm::Id > &connectivity, const std::string &coordsNm)
Definition: DataSetBuilderExplicit.h:312
vtkm::cont::ConvertNumComponentsToOffsets
void ConvertNumComponentsToOffsets(const vtkm::cont::UnknownArrayHandle &numComponentsArray, vtkm::cont::ArrayHandle< vtkm::Id > &offsetsArray, vtkm::Id &componentsArraySize, vtkm::cont::DeviceAdapterId device=vtkm::cont::DeviceAdapterTagAny{})
ConvertNumComponentsToOffsets takes an array of Vec sizes (i.e.
vtkm::cont::DataSetBuilderExplicitIterative::AddPoint
vtkm::Id AddPoint(const vtkm::Vec< T, 3 > &pt)
Add a point to the DataSet.
Definition: DataSetBuilderExplicit.h:393
vtkm::cont::DataSetBuilderExplicitIterative::connectivity
std::vector< vtkm::Id > connectivity
Definition: DataSetBuilderExplicit.h:464
VTKM_CONT_EXPORT
#define VTKM_CONT_EXPORT
Definition: vtkm_cont_export.h:44
VTKM_CONT
#define VTKM_CONT
Definition: ExportMacros.h:57
vtkm::Id
vtkm::Int64 Id
Base type to use to index arrays.
Definition: Types.h:227
vtkm::UInt8
uint8_t UInt8
Base type to use for 8-bit unsigned integer numbers.
Definition: Types.h:169
vtkm::CopyFlag::On
@ On
vtkm::cont::DataSet::AddCoordinateSystem
vtkm::IdComponent AddCoordinateSystem(const vtkm::cont::CoordinateSystem &cs)
Adds the given CoordinateSystem to the DataSet.
vtkm::cont::DataSetBuilderExplicit::DataSetBuilderExplicit
DataSetBuilderExplicit()
Definition: DataSetBuilderExplicit.h:31
vtkm::Vec< T, 3 >
Definition: Types.h:1016
vtkm::Vec< vtkm::FloatDefault, 3 >
vtkm::CopyFlag::Off
@ Off
vtkm::cont::DataSetBuilderExplicitIterative::numIdx
std::vector< vtkm::IdComponent > numIdx
Definition: DataSetBuilderExplicit.h:463
vtkm::FloatDefault
vtkm::Float32 FloatDefault
The floating point type to use when no other precision is specified.
Definition: Types.h:236
ConvertNumComponentsToOffsets.h
vtkm::cont::CellSetExplicit
Defines an irregular collection of cells.
Definition: CastAndCall.h:36
vtkm::cont::DataSetBuilderExplicitIterative::shapes
std::vector< vtkm::UInt8 > shapes
Definition: DataSetBuilderExplicit.h:462
vtkm::cont::ArrayHandle::Allocate
void Allocate(vtkm::Id numberOfValues, vtkm::CopyFlag preserve, vtkm::cont::Token &token) const
Allocates an array large enough to hold the given number of values.
Definition: ArrayHandle.h:490
vtkm::cont::DataSetBuilderExplicit
Definition: DataSetBuilderExplicit.h:27
vtkm::cont::make_ArrayHandle
vtkm::cont::ArrayHandleBasic< T > make_ArrayHandle(const T *array, vtkm::Id numberOfValues, vtkm::CopyFlag copy)
A convenience function for creating an ArrayHandle from a standard C array.
Definition: ArrayHandleBasic.h:270
vtkm::cont::ArrayHandle::WritePortal
WritePortalType WritePortal() const
Get an array portal that can be used in the control environment.
Definition: ArrayHandle.h:454
vtkm::cont::DataSetBuilderExplicitIterative::coordNm
std::string coordNm
Definition: DataSetBuilderExplicit.h:459
vtkm::cont::DataSetBuilderExplicit::Create
static vtkm::cont::DataSet Create(const std::vector< T > &xVals, const std::vector< T > &yVals, const std::vector< vtkm::UInt8 > &shapes, const std::vector< vtkm::IdComponent > &numIndices, const std::vector< vtkm::Id > &connectivity, const std::string &coordsNm="coords")
Create a 2D DataSet with arbitrary cell connectivity.
Definition: DataSetBuilderExplicit.h:83
vtkm::cont::DataSetBuilderExplicit::Create
static vtkm::cont::DataSet Create(const std::vector< T > &xVals, const std::vector< vtkm::UInt8 > &shapes, const std::vector< vtkm::IdComponent > &numIndices, const std::vector< vtkm::Id > &connectivity, const std::string &coordsNm="coords")
Create a 1D DataSet with arbitrary cell connectivity.
Definition: DataSetBuilderExplicit.h:52
DataSet.h
vtkm::cont::DataSetBuilderExplicitIterative
Helper class to build a DataSet by iteratively adding points and cells.
Definition: DataSetBuilderExplicit.h:369