VTK
/Users/kitware/Dashboards/MyTests/VTK_BLD_Release_docs/Utilities/Doxygen/dox/Accelerators/Dax/vtkDaxThresholdImpl.h
Go to the documentation of this file.
00001 //=============================================================================
00002 //
00003 //  Copyright (c) Kitware, Inc.
00004 //  All rights reserved.
00005 //  See LICENSE.txt for details.
00006 //
00007 //  This software is distributed WITHOUT ANY WARRANTY; without even
00008 //  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
00009 //  PURPOSE.  See the above copyright notice for more information.
00010 //
00011 //  Copyright 2012 Sandia Corporation.
00012 //  Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00013 //  the U.S. Government retains certain rights in this software.
00014 //
00015 //=============================================================================
00016 
00017 #ifndef vtkDaxThresholdImpl_h
00018 #define vtkDaxThresholdImpl_h
00019 
00020 // Common code
00021 #include "vtkDaxConfig.h"
00022 #include "vtkDaxDetailCommon.h"
00023 
00024 #include "vtkDispatcher.h"
00025 #include "vtkDoubleDispatcher.h"
00026 #include "vtkNew.h"
00027 
00028 //cell types we support
00029 #include "vtkCellTypes.h"
00030 #include "vtkGenericCell.h"
00031 #include "vtkHexahedron.h"
00032 #include "vtkLine.h"
00033 #include "vtkQuad.h"
00034 #include "vtkTetra.h"
00035 #include "vtkTriangle.h"
00036 #include "vtkVertex.h"
00037 #include "vtkVoxel.h"
00038 #include "vtkWedge.h"
00039 
00040 //fields we support
00041 #include "vtkDoubleArray.h"
00042 #include "vtkFloatArray.h"
00043 #include "vtkIntArray.h"
00044 #include "vtkUnsignedCharArray.h"
00045 
00046 //datasets we support
00047 #include "vtkDataObjectTypes.h"
00048 #include "vtkImageData.h"
00049 #include "vtkStructuredGrid.h"
00050 #include "vtkUniformGrid.h"
00051 #include "vtkUnstructuredGrid.h"
00052 
00053 //helpers that convert to and from Dax
00054 #include "vtkToDax/CellTypeToType.h"
00055 #include "vtkToDax/Containers.h"
00056 #include "vtkToDax/DataSetTypeToType.h"
00057 #include "vtkToDax/FieldTypeToType.h"
00058 #include "vtkToDax/Portals.h"
00059 #include "vtkToDax/Threshold.h"
00060 
00061 
00062 namespace vtkDax{
00063 namespace detail{
00064 
00065   struct ValidThresholdInput
00066   {
00067     typedef int ReturnType;
00068     vtkDataSet* Input;
00069     vtkCell* Cell;
00070     double Min;
00071     double Max;
00072 
00073     vtkUnstructuredGrid* Result;
00074 
00075     ValidThresholdInput(vtkDataSet* in, vtkUnstructuredGrid* out,
00076                         vtkCell* cell, double lower, double upper):
00077       Input(in),Cell(cell),Min(lower),Max(upper),Result(out){}
00078 
00079     template<typename LHS>
00080     int operator()(LHS &arrayField) const
00081       {
00082       //we can derive the type of the field at compile time, but not the
00083       //length
00084       switch(arrayField.GetNumberOfComponents())
00085         {
00086           case 1:
00087             //first we extract the field type of the array
00088             //second we extract the number of components
00089             typedef typename vtkToDax::FieldTypeToType<LHS,1>::DaxValueType VT1;
00090             return dispatchOnFieldType<LHS,VT1>(arrayField);
00091           case 2:
00092             typedef typename vtkToDax::FieldTypeToType<LHS,2>::DaxValueType VT2;
00093           return dispatchOnFieldType<LHS,VT2>(arrayField);
00094           case 3:
00095             typedef typename vtkToDax::FieldTypeToType<LHS,3>::DaxValueType VT3;
00096             return dispatchOnFieldType<LHS,VT3>(arrayField);
00097         default:
00098           //currently only support 1 to 3 components
00099           //we need to make dispatch on field data smarter in that it does
00100           //this automagically
00101           return 0;
00102         }
00103 
00104 
00105       }
00106 
00107     template<typename VTKArrayType, typename DaxValueType>
00108     int dispatchOnFieldType(VTKArrayType& vtkField) const
00109       {
00110       typedef vtkToDax::vtkArrayContainerTag<VTKArrayType> FieldTag;
00111       typedef dax::cont::ArrayHandle<DaxValueType,FieldTag> FieldHandle;
00112       typedef typename dax::cont::ArrayHandle<DaxValueType,
00113                       FieldTag>::PortalConstControl      PortalType;
00114 
00115       FieldHandle field = FieldHandle( PortalType(&vtkField,
00116                                             vtkField.GetNumberOfTuples() ) );
00117       vtkToDax::Threshold<FieldHandle> threshold(field,
00118                                                  DaxValueType(Min),
00119                                                  DaxValueType(Max));
00120       threshold.setFieldName(vtkField.GetName());
00121       threshold.setOutputGrid(this->Result);
00122 
00123       //see if we have a valid data set type
00124       //if so will perform the threshold if possible
00125       vtkDoubleDispatcher<vtkDataSet,vtkCell,int> dataDispatcher;
00126       dataDispatcher.Add<vtkImageData,vtkVoxel>(threshold);
00127       dataDispatcher.Add<vtkUniformGrid,vtkVoxel>(threshold);
00128 
00129       dataDispatcher.Add<vtkUnstructuredGrid,vtkHexahedron>(threshold);
00130       dataDispatcher.Add<vtkUnstructuredGrid,vtkLine>(threshold);
00131       dataDispatcher.Add<vtkUnstructuredGrid,vtkQuad>(threshold);
00132       dataDispatcher.Add<vtkUnstructuredGrid,vtkTetra>(threshold);
00133       dataDispatcher.Add<vtkUnstructuredGrid,vtkTriangle>(threshold);
00134       dataDispatcher.Add<vtkUnstructuredGrid,vtkVertex>(threshold);
00135       dataDispatcher.Add<vtkUnstructuredGrid,vtkWedge>(threshold);
00136 
00137       int validThreshold = dataDispatcher.Go(this->Input,this->Cell);
00138       return validThreshold;
00139       }
00140   private:
00141     void operator=(const ValidThresholdInput&);
00142   };
00143 } //end detail namespace
00144 
00145 
00146 //------------------------------------------------------------------------------
00147 int Threshold(vtkDataSet* input, vtkUnstructuredGrid *output,
00148               vtkDataArray* field, double lower, double upper)
00149 {
00150   //we are doing a point threshold now verify we have suitable cells
00151   //Dax currently supports: hexs,lines,quads,tets,triangles,vertex,voxel,wedge
00152   //if something a cell that doesn't match that list we punt to the
00153   //VTK implementation.
00154   vtkDax::detail::CellTypeInDataSet cType = vtkDax::detail::cellType(input);
00155 
00156   //construct the object that holds all the state needed to do the threshold
00157   vtkDax::detail::ValidThresholdInput validInput(input,output,cType.Cell,
00158                                                  lower,
00159                                                  upper);
00160 
00161 
00162   //setup the dispatch to only allow float and int array to go to the next step
00163   vtkDispatcher<vtkAbstractArray,int> fieldDispatcher;
00164   fieldDispatcher.Add<vtkFloatArray>(validInput);
00165   fieldDispatcher.Add<vtkDoubleArray>(validInput);
00166   fieldDispatcher.Add<vtkUnsignedCharArray>(validInput);
00167   fieldDispatcher.Add<vtkIntArray>(validInput);
00168   return fieldDispatcher.Go(field);
00169 }
00170 
00171 } //end vtkDax namespace
00172 // VTK-HeaderTest-Exclude: vtkDaxThresholdImpl.h
00173 #endif