VTK
|
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 vtkDaxContourImpl_h 00018 #define vtkDaxContourImpl_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 //fields we support 00029 #include "vtkDoubleArray.h" 00030 #include "vtkFloatArray.h" 00031 00032 //cell types we support 00033 #include "vtkCellTypes.h" 00034 #include "vtkGenericCell.h" 00035 #include "vtkHexahedron.h" 00036 #include "vtkTetra.h" 00037 #include "vtkTriangle.h" 00038 #include "vtkVoxel.h" 00039 00040 //datasets we support 00041 #include "vtkDataObjectTypes.h" 00042 #include "vtkImageData.h" 00043 #include "vtkUniformGrid.h" 00044 #include "vtkPolyData.h" 00045 00046 //helpers that convert vtk to dax 00047 #include "vtkToDax/CellTypeToType.h" 00048 #include "vtkToDax/Containers.h" 00049 #include "vtkToDax/Contour.h" 00050 #include "vtkToDax/DataSetTypeToType.h" 00051 #include "vtkToDax/FieldTypeToType.h" 00052 #include "vtkToDax/Portals.h" 00053 00054 namespace vtkDax { 00055 namespace detail { 00056 struct ValidContourInput 00057 { 00058 typedef int ReturnType; 00059 vtkDataSet* Input; 00060 vtkCell* Cell; 00061 double IsoValue; 00062 bool ComputeScalars; 00063 00064 vtkPolyData* Result; 00065 00066 ValidContourInput(vtkDataSet* in, vtkPolyData* out, 00067 vtkCell* cell, double isoValue, 00068 bool computeScalars) : 00069 Input(in), 00070 Cell(cell), 00071 IsoValue(isoValue), 00072 ComputeScalars(computeScalars), 00073 Result(out) { } 00074 00075 template<typename LHS> 00076 int operator()(LHS &arrayField) const 00077 { 00078 //we can derive the type of the field at compile time, but not the 00079 //length 00080 if (arrayField.GetNumberOfComponents() == 1) 00081 { 00082 //first we extract the field type of the array 00083 //second we extract the number of components 00084 typedef typename vtkToDax::FieldTypeToType<LHS,1>::DaxValueType VT1; 00085 return dispatchOnFieldType<LHS,VT1>(arrayField); 00086 } 00087 return 0; 00088 } 00089 00090 template<typename VTKArrayType, typename DaxValueType> 00091 int dispatchOnFieldType(VTKArrayType& vtkField) const 00092 { 00093 typedef vtkToDax::vtkArrayContainerTag<VTKArrayType> FieldTag; 00094 typedef dax::cont::ArrayHandle<DaxValueType,FieldTag> FieldHandle; 00095 00096 typedef typename dax::cont::ArrayHandle 00097 <DaxValueType, FieldTag>::PortalConstControl PortalType; 00098 00099 FieldHandle field = FieldHandle( PortalType(&vtkField, 00100 vtkField.GetNumberOfTuples() ) ); 00101 vtkToDax::Contour<FieldHandle> contour(field, 00102 DaxValueType(this->IsoValue), 00103 this->ComputeScalars); 00104 contour.setFieldName(vtkField.GetName()); 00105 contour.setOutputGrid(this->Result); 00106 00107 // see if we have a valid data set type if so will perform the 00108 // marchingcubes if possible 00109 vtkDoubleDispatcher<vtkDataSet,vtkCell,int> dataDispatcher; 00110 dataDispatcher.Add<vtkImageData,vtkVoxel>(contour); 00111 dataDispatcher.Add<vtkUniformGrid,vtkVoxel>(contour); 00112 dataDispatcher.Add<vtkUnstructuredGrid,vtkHexahedron>(contour); 00113 dataDispatcher.Add<vtkUnstructuredGrid,vtkTetra>(contour); 00114 00115 int validMC = dataDispatcher.Go(this->Input,this->Cell); 00116 return validMC; 00117 } 00118 private: 00119 void operator=(const ValidContourInput&); 00120 }; 00121 } //namespace detail 00122 00123 00124 //------------------------------------------------------------------------------ 00125 int Contour(vtkDataSet* input, vtkPolyData *output, 00126 vtkDataArray* field, float isoValue, bool computeScalars) 00127 { 00128 //we are doing a point threshold now verify we have suitable cells 00129 //Dax currently supports: hexs,lines,quads,tets,triangles,vertex,voxel,wedge 00130 //if something a cell that doesn't match that list we punt to the 00131 //VTK implementation. 00132 vtkDax::detail::CellTypeInDataSet cType = vtkDax::detail::cellType(input); 00133 00134 //construct the object that holds all the state needed to do the MC 00135 vtkDax::detail::ValidContourInput validInput(input,output,cType.Cell, 00136 isoValue, computeScalars); 00137 00138 00139 //setup the dispatch to only allow float and int array to go to the next step 00140 vtkDispatcher<vtkAbstractArray,int> fieldDispatcher; 00141 fieldDispatcher.Add<vtkFloatArray>(validInput); 00142 fieldDispatcher.Add<vtkDoubleArray>(validInput); 00143 return fieldDispatcher.Go(field); 00144 } 00145 00146 } //end vtkDax namespace 00147 // VTK-HeaderTest-Exclude: vtkDaxContourImpl.h 00148 #endif