![]() |
Mesh Oriented datABase
(version 5.4.1)
Array-based unstructured mesh datastructure
|
00001 #ifndef __smoab_DataSetConverter_h
00002 #define __smoab_DataSetConverter_h
00003
00004 #include "SimpleMoab.h"
00005 #include "CellSets.h"
00006 #include "detail/LoadGeometry.h"
00007 #include "detail/ReadSparseTag.h"
00008
00009 #include
00010 #include
00011 #include
00012 #include
00013 #include
00014 #include
00015
00016 #include
00017
00018 namespace smoab
00019 {
00020 class DataSetConverter
00021 {
00022 const smoab::Interface& Interface;
00023 moab::Interface* Moab;
00024 const smoab::Tag* Tag;
00025 bool ReadMaterialIds;
00026 bool ReadProperties;
00027 std::string MaterialName;
00028
00029 public:
00030 DataSetConverter( const smoab::Interface& interface, const smoab::Tag* tag )
00031 : Interface( interface ), Moab( interface.Moab ), Tag( tag ), ReadMaterialIds( false ), ReadProperties( false )
00032 {
00033 }
00034
00035 void readMaterialIds( bool add )
00036 {
00037 this->ReadMaterialIds = add;
00038 }
00039 bool readMaterialIds() const
00040 {
00041 return this->ReadMaterialIds;
00042 }
00043
00044 void readProperties( bool readProps )
00045 {
00046 this->ReadProperties = readProps;
00047 }
00048 bool readProperties() const
00049 {
00050 return this->ReadProperties;
00051 }
00052
00053 //----------------------------------------------------------------------------
00054 //given a range of entity handles merge them all into a single unstructured
00055 //grid. Currently doesn't support reading properties.
00056 //Will read in material ids, if no material id is assigned to an entity,
00057 //its cells will be given an unique id
00058 template < typename VTKGridType >
00059 bool fill( const smoab::Range& entities, VTKGridType* grid ) const
00060 {
00061 //create a helper datastructure which can determines all the unique point ids
00062 //and converts moab connecitvity info to vtk connectivity
00063
00064 //get all the cells for each parent entity and create
00065 // an entity set of those items
00066 int dim = this->Tag->value();
00067 typedef smoab::Range::const_iterator iterator;
00068 smoab::CellSets entitySets;
00069 for( iterator i = entities.begin(); i != entities.end(); ++i )
00070 {
00071 smoab::Range entitiesCells;
00072 if( this->Tag->isComparable() )
00073 {
00074 //if we are comparable only find the cells that match our tags dimension
00075 entitiesCells = this->Interface.findEntitiesWithDimension( *i, dim, true );
00076 }
00077 else
00078 {
00079 entitiesCells = this->Interface.findHighestDimensionEntities( *i, true );
00080 }
00081 smoab::CellSet set( *i, entitiesCells );
00082 entitySets.push_back( set );
00083 }
00084
00085 moab::Range cells = smoab::getAllCells( entitySets );
00086
00087 //convert the datastructure from a list of cells to a vtk data set
00088 detail::LoadGeometry loadGeom( cells, dim, this->Interface );
00089 loadGeom.fill( grid );
00090
00091 if( this->readMaterialIds() )
00092 {
00093
00094 detail::ReadSparseTag materialTagReading( entitySets, cells, this->Interface );
00095
00096 smoab::MaterialTag mtag;
00097 vtkNew< vtkIntArray > materials;
00098 materials->SetName( mtag.name() );
00099 materialTagReading.fill( materials.GetPointer(), &mtag );
00100 grid->GetCellData()->AddArray( materials.GetPointer() );
00101 }
00102
00103 //by default we always try to load the default tag
00104 detail::ReadSparseTag sTagReading( entitySets, cells, this->Interface );
00105
00106 vtkNew< vtkIntArray > sparseTagData;
00107 sparseTagData->SetName( this->Tag->name() );
00108 sTagReading.fill( sparseTagData.GetPointer(), this->Tag );
00109 grid->GetCellData()->AddArray( sparseTagData.GetPointer() );
00110
00111 return true;
00112 }
00113
00114 //----------------------------------------------------------------------------
00115 //given a single entity handle create a unstructured grid from it.
00116 //optional third parameter is the material id to use if readMaterialIds
00117 //is on, and no material sparse tag is found for this entity
00118 template < typename VTKGridType >
00119 bool fill( const smoab::EntityHandle& entity, VTKGridType* grid, const int materialId = 0 ) const
00120 {
00121 //create a helper datastructure which can determines all the unique point ids
00122 //and converts moab connecitvity info to vtk connectivity
00123
00124 smoab::Range cells;
00125 int dim = this->Tag->value();
00126 if( this->Tag->isComparable() )
00127 {
00128 //if we are comparable only find the cells that match our tags dimension
00129 cells = this->Interface.findEntitiesWithDimension( entity, dim, true );
00130 }
00131 else
00132 {
00133 //load subentities
00134 cells = this->Interface.findHighestDimensionEntities( entity, true );
00135 }
00136
00137 //convert the datastructure from a list of cells to a vtk data set
00138 detail::LoadGeometry loadGeom( cells, dim, this->Interface );
00139 loadGeom.fill( grid );
00140
00141 const smoab::Range& points = loadGeom.moabPoints();
00142
00143 if( this->readProperties() )
00144 {
00145 this->readProperties( cells, grid->GetCellData() );
00146 this->readProperties( points, grid->GetPointData() );
00147 }
00148
00149 smoab::CellSets cellSets;
00150 smoab::CellSet set( entity, cells );
00151 cellSets.push_back( set );
00152 if( this->readMaterialIds() )
00153 {
00154 smoab::MaterialTag mtag;
00155 detail::ReadSparseTag materialTagReading( cellSets, cells, this->Interface );
00156
00157 vtkNew< vtkIntArray > materials;
00158 materials->SetName( mtag.name() );
00159 materialTagReading.fill( materials.GetPointer(), &mtag );
00160 grid->GetCellData()->AddArray( materials.GetPointer() );
00161 }
00162
00163 //by default we always try to load the default tag
00164 detail::ReadSparseTag sTagReading( cellSets, cells, this->Interface );
00165
00166 vtkNew< vtkIntArray > sparseTagData;
00167 sparseTagData->SetName( this->Tag->name() );
00168 sTagReading.fill( sparseTagData.GetPointer(), this->Tag );
00169 grid->GetCellData()->AddArray( sparseTagData.GetPointer() );
00170
00171 return true;
00172 }
00173
00174 private:
00175 //----------------------------------------------------------------------------
00176 void readProperties( smoab::Range const& entities, vtkFieldData* field ) const
00177 {
00178 if( entities.empty() )
00179 {
00180 return;
00181 }
00182
00183 //so we get all the tags and parse out the sparse and dense tags
00184 //that we support
00185 typedef std::vector< moab::Tag >::const_iterator iterator;
00186 std::vector< moab::Tag > tags;
00187 this->Moab->tag_get_tags_on_entity( entities.front(), tags );
00188
00189 this->readDenseTags( tags, entities, field );
00190 }
00191
00192 //----------------------------------------------------------------------------
00193 void readDenseTags( std::vector< moab::Tag >& tags, smoab::Range const& entities, vtkFieldData* field ) const
00194 {
00195 typedef std::vector< moab::Tag >::const_iterator iterator;
00196
00197 for( iterator i = tags.begin(); i != tags.end(); ++i )
00198 {
00199 moab::TagType tagType;
00200 moab::DataType tagDataType;
00201
00202 this->Moab->tag_get_type( *i, tagType );
00203 this->Moab->tag_get_data_type( *i, tagDataType );
00204
00205 //make sure it is only dense
00206 if( tagType != moab::MB_TAG_DENSE )
00207 {
00208 continue;
00209 }
00210 //and only integer and double
00211 if( tagDataType != moab::MB_TYPE_DOUBLE && tagDataType != moab::MB_TYPE_INTEGER )
00212 {
00213 //unsupported type, skip to next tag
00214 continue;
00215 }
00216
00217 //read the name of the tag
00218 std::string name;
00219 name.reserve( 32 );
00220 this->Moab->tag_get_name( *i, name );
00221
00222 //read the number of components of the tag
00223 int numComps = 1;
00224
00225 this->Moab->tag_get_length( *i, numComps );
00226
00227 //read the data if it is one of the two types we support
00228 int size = entities.size();
00229 if( tagDataType == moab::MB_TYPE_DOUBLE )
00230 {
00231 vtkNew< vtkDoubleArray > array;
00232 array->SetName( name.c_str() );
00233 array->SetNumberOfComponents( numComps );
00234 array->SetNumberOfTuples( size );
00235
00236 //read directly into the double array
00237 this->Moab->tag_get_data( *i, entities, array->GetVoidPointer( 0 ) );
00238 field->AddArray( array.GetPointer() );
00239 }
00240 else if( tagDataType == moab::MB_TYPE_INTEGER )
00241 {
00242 vtkNew< vtkIntArray > array;
00243 array->SetName( name.c_str() );
00244 array->SetNumberOfComponents( numComps );
00245 array->SetNumberOfTuples( size );
00246
00247 //read directly into the double array
00248 this->Moab->tag_get_data( *i, entities, array->GetVoidPointer( 0 ) );
00249 field->AddArray( array.GetPointer() );
00250 }
00251 else
00252 {
00253 }
00254 }
00255 }
00256 };
00257
00258 } // namespace smoab
00259 #endif // __smoab_DataSetConverter_h