MOAB: Mesh Oriented datABase  (version 5.4.1)
LoadGeometry.h
Go to the documentation of this file.
00001 #ifndef __smoab_detail_LoadGeometry_h
00002 #define __smoab_detail_LoadGeometry_h
00003 
00004 #include "SimpleMoab.h"
00005 #include "MixedCellConnectivity.h"
00006 
00007 #include <vtkCellArray.h>
00008 #include <vtkIdTypeArray.h>
00009 #include <vtkNew.h>
00010 #include <vtkPointData.h>
00011 #include <vtkPoints.h>
00012 #include <vtkPointSet.h>
00013 #include <vtkPolyData.h>
00014 #include <vtkUnsignedCharArray.h>
00015 #include <vtkUnstructuredGrid.h>
00016 namespace smoab
00017 {
00018 namespace detail
00019 {
00020 
00021     namespace
00022     {
00023 
00024         //empty definition, only allow specializations to compile
00025         template < typename VTKGridType >
00026         struct setCells;
00027 
00028         template <>
00029         struct setCells< vtkUnstructuredGrid >
00030         {
00031             typedef vtkUnstructuredGrid VTKGridType;
00032             void operator()( VTKGridType* grid,
00033                              int,
00034                              vtkUnsignedCharArray* types,
00035                              vtkIdTypeArray* locations,
00036                              vtkCellArray* cells ) const
00037             {
00038                 grid->SetCells( types, locations, cells, NULL, NULL );
00039             }
00040         };
00041 
00042         template <>
00043         struct setCells< vtkPolyData >
00044         {
00045             typedef vtkPolyData VTKGridType;
00046             void operator()( VTKGridType* grid,
00047                              int dim,
00048                              vtkUnsignedCharArray*,
00049                              vtkIdTypeArray*,
00050                              vtkCellArray* cells ) const
00051             {
00052                 if( dim == 0 )
00053                 {
00054                     grid->SetVerts( cells );
00055                 }
00056                 else if( dim == 1 )
00057                 {
00058                     grid->SetLines( cells );
00059                 }
00060                 else if( dim == 2 )
00061                 {
00062                     grid->SetPolys( cells );
00063                 }
00064             }
00065         };
00066 
00067     }  // namespace
00068     //basic class that given a range of moab cells will convert them to am unstructured grid.
00069     //holds only references to input cells, so they can't go out of scope
00070 
00071     //The Topology tag is used to describe the dimensoniality of the cells we are reading
00072     // 0 == verts, 1 == lines, 2 == tri/quads, 3 = volume elements
00073     class LoadGeometry
00074     {
00075         moab::Interface* Interface;
00076         int TopologyDim;
00077         smoab::detail::MixedCellConnectivity MixConn;
00078         smoab::Range Points;
00079 
00080       public:
00081         //----------------------------------------------------------------------------
00082         //warning cells input to constructor is held by reference
00083         LoadGeometry( const smoab::Range& cells, int topologyDim, const smoab::Interface& interface )
00084             : Interface( interface.Moab ), TopologyDim( topologyDim ), MixConn( cells, interface.Moab ), Points()
00085         {
00086         }
00087 
00088         const smoab::Range& moabPoints() const
00089         {
00090             return this->Points;
00091         }
00092 
00093         //----------------------------------------------------------------------------
00094         //todo: have support for using only a subsection of the input cells
00095         template < typename vtkDataSetType >
00096         void fill( vtkDataSetType* dataSet )
00097         {
00098             //now that mixConn has all the cells properly stored, lets fixup
00099             //the ids so that they start at zero and keep the same logical ordering
00100             //as before.
00101             vtkIdType numCells, connLen;
00102             this->MixConn.compactIds( numCells, connLen );
00103             this->addGridsTopology( dataSet, numCells, connLen );
00104             this->addCoordinates( dataSet );
00105         }
00106 
00107       private:
00108         //----------------------------------------------------------------------------
00109         void addCoordinates( vtkPointSet* grid )
00110         {
00111             //this is sorta of hackish as moabPoints is only valid
00112             //after compactIds has been called
00113             this->MixConn.moabPoints( this->Points );
00114 
00115             //since the smoab::range are always unique and sorted
00116             //we can use the more efficient coords_iterate
00117             //call in moab, which returns moab internal allocated memory
00118             vtkNew< vtkPoints > newPoints;
00119             newPoints->SetDataTypeToDouble();
00120             newPoints->SetNumberOfPoints( this->Points.size() );
00121 
00122             //need a pointer to the allocated vtkPoints memory so that we
00123             //don't need to use an extra copy and we can bypass all vtk's check
00124             //on out of bounds
00125             double* rawPoints = static_cast< double* >( newPoints->GetVoidPointer( 0 ) );
00126             this->Interface->get_coords( this->Points, rawPoints );
00127 
00128             grid->SetPoints( newPoints.GetPointer() );
00129         }
00130 
00131         //----------------------------------------------------------------------------
00132         template < typename vtkDataSetType >
00133         void addGridsTopology( vtkDataSetType* grid, vtkIdType numCells, vtkIdType numConnectivity ) const
00134         {
00135             //correct the connectivity size to account for the vtk padding
00136             const vtkIdType vtkConnectivity = numCells + numConnectivity;
00137 
00138             vtkNew< vtkIdTypeArray > cellArray;
00139             vtkNew< vtkIdTypeArray > cellLocations;
00140             vtkNew< vtkUnsignedCharArray > cellTypes;
00141 
00142             cellArray->SetNumberOfValues( vtkConnectivity );
00143             cellLocations->SetNumberOfValues( numCells );
00144             cellTypes->SetNumberOfValues( numCells );
00145 
00146             vtkIdType* rawArray     = static_cast< vtkIdType* >( cellArray->GetVoidPointer( 0 ) );
00147             vtkIdType* rawLocations = static_cast< vtkIdType* >( cellLocations->GetVoidPointer( 0 ) );
00148             unsigned char* rawTypes = static_cast< unsigned char* >( cellTypes->GetVoidPointer( 0 ) );
00149 
00150             this->MixConn.copyToVtkCellInfo( rawArray, rawLocations, rawTypes );
00151 
00152             vtkNew< vtkCellArray > cells;
00153             cells->SetCells( numCells, cellArray.GetPointer() );
00154 
00155             setCells< vtkDataSetType >()( grid, this->TopologyDim, cellTypes.GetPointer(), cellLocations.GetPointer(),
00156                                           cells.GetPointer() );
00157         }
00158     };
00159 
00160 }  // namespace detail
00161 }  // namespace smoab
00162 
00163 #endif  // __smoab_detail_LoadGeometry_h
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines