VTK
vtkDaxThresholdImpl.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 vtkDaxThresholdImpl_h
18 #define vtkDaxThresholdImpl_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 //cell types we support
29 #include "vtkCellTypes.h"
30 #include "vtkGenericCell.h"
31 #include "vtkHexahedron.h"
32 #include "vtkLine.h"
33 #include "vtkQuad.h"
34 #include "vtkTetra.h"
35 #include "vtkTriangle.h"
36 #include "vtkVertex.h"
37 #include "vtkVoxel.h"
38 #include "vtkWedge.h"
39 
40 //fields we support
41 #include "vtkDoubleArray.h"
42 #include "vtkFloatArray.h"
43 #include "vtkIntArray.h"
44 #include "vtkUnsignedCharArray.h"
45 
46 //datasets we support
47 #include "vtkDataObjectTypes.h"
48 #include "vtkImageData.h"
49 #include "vtkStructuredGrid.h"
50 #include "vtkUniformGrid.h"
51 #include "vtkUnstructuredGrid.h"
52 
53 //helpers that convert to and from Dax
55 #include "vtkToDax/Containers.h"
58 #include "vtkToDax/Portals.h"
59 #include "vtkToDax/Threshold.h"
60 
61 
62 namespace vtkDax{
63 namespace detail{
64 
66  {
67  typedef int ReturnType;
70  double Min;
71  double Max;
72 
74 
76  vtkCell* cell, double lower, double upper):
77  Input(in),Cell(cell),Min(lower),Max(upper),Result(out){}
78 
79  template<typename LHS>
80  int operator()(LHS &arrayField) const
81  {
82  //we can derive the type of the field at compile time, but not the
83  //length
84  switch(arrayField.GetNumberOfComponents())
85  {
86  case 1:
87  //first we extract the field type of the array
88  //second we extract the number of components
90  return dispatchOnFieldType<LHS,VT1>(arrayField);
91  case 2:
93  return dispatchOnFieldType<LHS,VT2>(arrayField);
94  case 3:
96  return dispatchOnFieldType<LHS,VT3>(arrayField);
97  default:
98  //currently only support 1 to 3 components
99  //we need to make dispatch on field data smarter in that it does
100  //this automagically
101  return 0;
102  }
103 
104 
105  }
106 
107  template<typename VTKArrayType, typename DaxValueType>
108  int dispatchOnFieldType(VTKArrayType& vtkField) const
109  {
111  typedef dax::cont::ArrayHandle<DaxValueType,FieldTag> FieldHandle;
112  typedef typename dax::cont::ArrayHandle<DaxValueType,
113  FieldTag>::PortalConstControl PortalType;
114 
115  FieldHandle field = FieldHandle( PortalType(&vtkField,
116  vtkField.GetNumberOfTuples() ) );
117  vtkToDax::Threshold<FieldHandle> threshold(field,
118  DaxValueType(Min),
119  DaxValueType(Max));
120  threshold.setFieldName(vtkField.GetName());
121  threshold.setOutputGrid(this->Result);
122 
123  //see if we have a valid data set type
124  //if so will perform the threshold if possible
126  dataDispatcher.Add<vtkImageData,vtkVoxel>(threshold);
127  dataDispatcher.Add<vtkUniformGrid,vtkVoxel>(threshold);
128 
129  dataDispatcher.Add<vtkUnstructuredGrid,vtkHexahedron>(threshold);
130  dataDispatcher.Add<vtkUnstructuredGrid,vtkLine>(threshold);
131  dataDispatcher.Add<vtkUnstructuredGrid,vtkQuad>(threshold);
132  dataDispatcher.Add<vtkUnstructuredGrid,vtkTetra>(threshold);
133  dataDispatcher.Add<vtkUnstructuredGrid,vtkTriangle>(threshold);
134  dataDispatcher.Add<vtkUnstructuredGrid,vtkVertex>(threshold);
135  dataDispatcher.Add<vtkUnstructuredGrid,vtkWedge>(threshold);
136 
137  int validThreshold = dataDispatcher.Go(this->Input,this->Cell);
138  return validThreshold;
139  }
140  private:
141  void operator=(const ValidThresholdInput&);
142  };
143 } //end detail namespace
144 
145 
146 //------------------------------------------------------------------------------
148  vtkDataArray* field, double lower, double upper)
149 {
150  //we are doing a point threshold now verify we have suitable cells
151  //Dax currently supports: hexs,lines,quads,tets,triangles,vertex,voxel,wedge
152  //if something a cell that doesn't match that list we punt to the
153  //VTK implementation.
155 
156  //construct the object that holds all the state needed to do the threshold
157  vtkDax::detail::ValidThresholdInput validInput(input,output,cType.Cell,
158  lower,
159  upper);
160 
161 
162  //setup the dispatch to only allow float and int array to go to the next step
163  vtkDispatcher<vtkAbstractArray,int> fieldDispatcher;
164  fieldDispatcher.Add<vtkFloatArray>(validInput);
165  fieldDispatcher.Add<vtkDoubleArray>(validInput);
166  fieldDispatcher.Add<vtkUnsignedCharArray>(validInput);
167  fieldDispatcher.Add<vtkIntArray>(validInput);
168  return fieldDispatcher.Go(field);
169 }
170 
171 } //end vtkDax namespace
172 // VTK-HeaderTest-Exclude: vtkDaxThresholdImpl.h
173 #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...
void Add(Functor fun)
Add in a functor that is mapped to the template SomeLhs parameter.
int operator()(LHS &arrayField) const
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
a cell that represents a 3D point
Definition: vtkVertex.h:36
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:41
void setFieldName(const char *name)
Definition: Threshold.h:180
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:41
int dispatchOnFieldType(VTKArrayType &vtkField) const
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
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:47
cell represents a 1D line
Definition: vtkLine.h:35
abstract class to specify cell behavior
Definition: vtkCell.h:59
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:45
Dispatch to functor based on two pointer types.
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:44
int Threshold(vtkDataSet *input, vtkUnstructuredGrid *output, vtkDataArray *field, double lower, double upper)
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
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:47
ValidThresholdInput(vtkDataSet *in, vtkUnstructuredGrid *out, vtkCell *cell, double lower, double upper)
Dispatch to functor based on a pointer type.
Definition: vtkDispatcher.h:91
dynamic, self-adjusting array of unsigned char
void setOutputGrid(vtkUnstructuredGrid *grid)
Definition: Threshold.h:175
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...
a cell that represents a triangle
Definition: vtkTriangle.h:41
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:49