VTK
vtkDaxContourImpl.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 vtkDaxContourImpl_h
18 #define vtkDaxContourImpl_h
19 
20 // Common code
21 #include "vtkDaxConfig.h"
22 #include "vtkDaxDetailCommon.h"
23 
24 #include "vtkDispatcher.h"
25 #include "vtkDoubleDispatcher.h"
26 #include "vtkNew.h"
27 
28 //fields we support
29 #include "vtkDoubleArray.h"
30 #include "vtkFloatArray.h"
31 
32 //cell types we support
33 #include "vtkCellTypes.h"
34 #include "vtkGenericCell.h"
35 #include "vtkHexahedron.h"
36 #include "vtkTetra.h"
37 #include "vtkTriangle.h"
38 #include "vtkVoxel.h"
39 
40 //datasets we support
41 #include "vtkDataObjectTypes.h"
42 #include "vtkImageData.h"
43 #include "vtkUniformGrid.h"
44 #include "vtkPolyData.h"
45 
46 //helpers that convert vtk to dax
48 #include "vtkToDax/Containers.h"
49 #include "vtkToDax/Contour.h"
52 #include "vtkToDax/Portals.h"
53 
54 namespace vtkDax {
55 namespace detail {
57  {
58  typedef int ReturnType;
61  double IsoValue;
63 
65 
67  vtkCell* cell, double isoValue,
68  bool computeScalars) :
69  Input(in),
70  Cell(cell),
71  IsoValue(isoValue),
72  ComputeScalars(computeScalars),
73  Result(out) { }
74 
75  template<typename LHS>
76  int operator()(LHS &arrayField) const
77  {
78  //we can derive the type of the field at compile time, but not the
79  //length
80  if (arrayField.GetNumberOfComponents() == 1)
81  {
82  //first we extract the field type of the array
83  //second we extract the number of components
85  return dispatchOnFieldType<LHS,VT1>(arrayField);
86  }
87  return 0;
88  }
89 
90  template<typename VTKArrayType, typename DaxValueType>
91  int dispatchOnFieldType(VTKArrayType& vtkField) const
92  {
94  typedef dax::cont::ArrayHandle<DaxValueType,FieldTag> FieldHandle;
95 
96  typedef typename dax::cont::ArrayHandle
97  <DaxValueType, FieldTag>::PortalConstControl PortalType;
98 
99  FieldHandle field = FieldHandle( PortalType(&vtkField,
100  vtkField.GetNumberOfTuples() ) );
101  vtkToDax::Contour<FieldHandle> contour(field,
102  DaxValueType(this->IsoValue),
103  this->ComputeScalars);
104  contour.setFieldName(vtkField.GetName());
105  contour.setOutputGrid(this->Result);
106 
107  // see if we have a valid data set type if so will perform the
108  // marchingcubes if possible
110  dataDispatcher.Add<vtkImageData,vtkVoxel>(contour);
111  dataDispatcher.Add<vtkUniformGrid,vtkVoxel>(contour);
112  dataDispatcher.Add<vtkUnstructuredGrid,vtkHexahedron>(contour);
113  dataDispatcher.Add<vtkUnstructuredGrid,vtkTetra>(contour);
114 
115  int validMC = dataDispatcher.Go(this->Input,this->Cell);
116  return validMC;
117  }
118  private:
119  void operator=(const ValidContourInput&);
120  };
121 } //namespace detail
122 
123 
124 //------------------------------------------------------------------------------
125 int Contour(vtkDataSet* input, vtkPolyData *output,
126  vtkDataArray* field, float isoValue, bool computeScalars)
127 {
128  //we are doing a point threshold now verify we have suitable cells
129  //Dax currently supports: hexs,lines,quads,tets,triangles,vertex,voxel,wedge
130  //if something a cell that doesn't match that list we punt to the
131  //VTK implementation.
133 
134  //construct the object that holds all the state needed to do the MC
135  vtkDax::detail::ValidContourInput validInput(input,output,cType.Cell,
136  isoValue, computeScalars);
137 
138 
139  //setup the dispatch to only allow float and int array to go to the next step
140  vtkDispatcher<vtkAbstractArray,int> fieldDispatcher;
141  fieldDispatcher.Add<vtkFloatArray>(validInput);
142  fieldDispatcher.Add<vtkDoubleArray>(validInput);
143  return fieldDispatcher.Go(field);
144 }
145 
146 } //end vtkDax namespace
147 // VTK-HeaderTest-Exclude: vtkDaxContourImpl.h
148 #endif
ReturnType Go(BaseLhs *lhs, BaseRhs *rhs)
Given two pointers of objects that derive from the BaseLhs and BaseRhs we find the matching functor t...
ValidContourInput(vtkDataSet *in, vtkPolyData *out, vtkCell *cell, double isoValue, bool computeScalars)
void Add(Functor fun)
Add in a functor that is mapped to the template SomeLhs parameter.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:41
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
void Add(Functor fun)
Add in a functor that is mapped to the combination of the two template parameters passed in...
dynamic, self-adjusting array of double
void setFieldName(const char *name)
Definition: Contour.h:246
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:47
abstract class to specify cell behavior
Definition: vtkCell.h:59
Dispatch to functor based on two pointer types.
int operator()(LHS &arrayField) const
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:44
topologically and geometrically regular array of data
Definition: vtkImageData.h:45
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
void setOutputGrid(vtkPolyData *grid)
Definition: Contour.h:241
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:47
Dispatch to functor based on a pointer type.
Definition: vtkDispatcher.h:91
CellTypeInDataSet cellType(vtkDataSet *input)
image data with blanking
ReturnType Go(BaseLhs *lhs)
Given a pointer to an object that derives from the BaseLhs we find the matching functor that was adde...
int Contour(vtkDataSet *input, vtkPolyData *output, vtkDataArray *field, float isoValue, bool computeScalars)
int dispatchOnFieldType(VTKArrayType &vtkField) const