Mesh Oriented datABase
(version 5.4.1)
Array-based unstructured mesh datastructure
|
00001 /** 00002 * MOAB, a Mesh-Oriented datABase, is a software component for creating, 00003 * storing and accessing finite element mesh data. 00004 * 00005 * Copyright 2004 Sandia Corporation. Under the terms of Contract 00006 * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government 00007 * retains certain rights in this software. 00008 * 00009 * This library is free software; you can redistribute it and/or 00010 * modify it under the terms of the GNU Lesser General Public 00011 * License as published by the Free Software Foundation; either 00012 * version 2.1 of the License, or (at your option) any later version. 00013 * 00014 */ 00015 00016 #ifdef WIN32 00017 #ifdef _DEBUG 00018 // turn off warnings that say they debugging identifier has been truncated 00019 // this warning comes up when using some STL containers 00020 #pragma warning( disable : 4786 ) 00021 #endif 00022 #endif 00023 00024 #include "WriteVtk.hpp" 00025 #include "moab/VtkUtil.hpp" 00026 #include "SysUtil.hpp" 00027 00028 #include <fstream> 00029 #include <iostream> 00030 #include <cstdio> 00031 #include <cassert> 00032 #include <vector> 00033 #include <set> 00034 #include <map> 00035 #include <iterator> 00036 00037 #include "moab/Interface.hpp" 00038 #include "moab/Range.hpp" 00039 #include "moab/CN.hpp" 00040 #include "MBTagConventions.hpp" 00041 #include "moab/WriteUtilIface.hpp" 00042 #include "Internals.hpp" 00043 #include "moab/FileOptions.hpp" 00044 00045 #define INS_ID( stringvar, prefix, id ) sprintf( stringvar, prefix, id ) 00046 00047 namespace moab 00048 { 00049 00050 const int DEFAULT_PRECISION = 10; 00051 const bool DEFAULT_STRICT = true; 00052 00053 WriterIface* WriteVtk::factory( Interface* iface ) 00054 { 00055 return new WriteVtk( iface ); 00056 } 00057 00058 WriteVtk::WriteVtk( Interface* impl ) 00059 : mbImpl( impl ), writeTool( 0 ), mStrict( DEFAULT_STRICT ), freeNodes( 0 ), createOneNodeCells( false ) 00060 { 00061 assert( impl != NULL ); 00062 impl->query_interface( writeTool ); 00063 } 00064 00065 WriteVtk::~WriteVtk() 00066 { 00067 mbImpl->release_interface( writeTool ); 00068 } 00069 00070 ErrorCode WriteVtk::write_file( const char* file_name, 00071 const bool overwrite, 00072 const FileOptions& opts, 00073 const EntityHandle* output_list, 00074 const int num_sets, 00075 const std::vector< std::string >& /* qa_list */, 00076 const Tag* tag_list, 00077 int num_tags, 00078 int /* export_dimension */ ) 00079 { 00080 ErrorCode rval; 00081 00082 // Get precision for node coordinates 00083 int precision; 00084 if( MB_SUCCESS != opts.get_int_option( "PRECISION", precision ) ) precision = DEFAULT_PRECISION; 00085 00086 if( MB_SUCCESS == opts.get_null_option( "STRICT" ) ) 00087 mStrict = true; 00088 else if( MB_SUCCESS == opts.get_null_option( "RELAXED" ) ) 00089 mStrict = false; 00090 else 00091 mStrict = DEFAULT_STRICT; 00092 00093 if( MB_SUCCESS == opts.get_null_option( "CREATE_ONE_NODE_CELLS" ) ) createOneNodeCells = true; 00094 00095 // Get entities to write 00096 Range nodes, elems; 00097 rval = gather_mesh( output_list, num_sets, nodes, elems ); 00098 if( MB_SUCCESS != rval ) return rval; 00099 00100 // Honor overwrite flag 00101 if( !overwrite ) 00102 { 00103 rval = writeTool->check_doesnt_exist( file_name ); 00104 if( MB_SUCCESS != rval ) return rval; 00105 } 00106 00107 // Create file 00108 std::ofstream file( file_name ); 00109 if( !file ) 00110 { 00111 MB_SET_ERR( MB_FILE_WRITE_ERROR, "Could not open file: " << file_name ); 00112 } 00113 file.precision( precision ); 00114 00115 // Write file 00116 if( ( rval = write_header( file ) ) != MB_SUCCESS || ( rval = write_nodes( file, nodes ) ) != MB_SUCCESS || 00117 ( rval = write_elems( file, nodes, elems ) ) != MB_SUCCESS || 00118 ( rval = write_tags( file, true, nodes, tag_list, num_tags ) ) != MB_SUCCESS || 00119 ( rval = write_tags( file, false, elems, tag_list, num_tags ) ) != MB_SUCCESS ) 00120 { 00121 file.close(); 00122 remove( file_name ); 00123 return rval; 00124 } 00125 00126 return MB_SUCCESS; 00127 } 00128 00129 ErrorCode WriteVtk::gather_mesh( const EntityHandle* set_list, int num_sets, Range& nodes, Range& elems ) 00130 { 00131 ErrorCode rval; 00132 int e; 00133 00134 if( !set_list || !num_sets ) 00135 { 00136 Range a; 00137 rval = mbImpl->get_entities_by_handle( 0, a ); 00138 if( MB_SUCCESS != rval ) return rval; 00139 00140 Range::const_iterator node_i, elem_i, set_i; 00141 node_i = a.lower_bound( a.begin(), a.end(), CREATE_HANDLE( MBVERTEX, 0, e ) ); 00142 elem_i = a.lower_bound( node_i, a.end(), CREATE_HANDLE( MBEDGE, 0, e ) ); 00143 set_i = a.lower_bound( elem_i, a.end(), CREATE_HANDLE( MBENTITYSET, 0, e ) ); 00144 nodes.merge( node_i, elem_i ); 00145 elems.merge( elem_i, set_i ); 00146 00147 // Filter out unsupported element types 00148 EntityType et = MBEDGE; 00149 for( et++; et < MBENTITYSET; et++ ) 00150 { 00151 if( VtkUtil::get_vtk_type( et, CN::VerticesPerEntity( et ) ) ) continue; 00152 Range::iterator eit = elems.lower_bound( elems.begin(), elems.end(), CREATE_HANDLE( et, 0, e ) ), 00153 ep1it = elems.lower_bound( elems.begin(), elems.end(), CREATE_HANDLE( et + 1, 0, e ) ); 00154 elems.erase( eit, ep1it ); 00155 } 00156 } 00157 else 00158 { 00159 std::set< EntityHandle > visited; 00160 std::vector< EntityHandle > sets; 00161 sets.reserve( num_sets ); 00162 std::copy( set_list, set_list + num_sets, std::back_inserter( sets ) ); 00163 while( !sets.empty() ) 00164 { 00165 // Get next set 00166 EntityHandle set = sets.back(); 00167 sets.pop_back(); 00168 // Skip sets we've already done 00169 if( !visited.insert( set ).second ) continue; 00170 00171 Range a; 00172 rval = mbImpl->get_entities_by_handle( set, a ); 00173 if( MB_SUCCESS != rval ) return rval; 00174 00175 Range::const_iterator node_i, elem_i, set_i; 00176 node_i = a.lower_bound( a.begin(), a.end(), CREATE_HANDLE( MBVERTEX, 0, e ) ); 00177 elem_i = a.lower_bound( node_i, a.end(), CREATE_HANDLE( MBEDGE, 0, e ) ); 00178 set_i = a.lower_bound( elem_i, a.end(), CREATE_HANDLE( MBENTITYSET, 0, e ) ); 00179 nodes.merge( node_i, elem_i ); 00180 elems.merge( elem_i, set_i ); 00181 std::copy( set_i, a.end(), std::back_inserter( sets ) ); 00182 00183 a.clear(); 00184 rval = mbImpl->get_child_meshsets( set, a ); 00185 std::copy( a.begin(), a.end(), std::back_inserter( sets ) ); 00186 } 00187 00188 for( Range::const_iterator ei = elems.begin(); ei != elems.end(); ++ei ) 00189 { 00190 std::vector< EntityHandle > connect; 00191 rval = mbImpl->get_connectivity( &( *ei ), 1, connect ); 00192 if( MB_SUCCESS != rval ) return rval; 00193 00194 for( unsigned int i = 0; i < connect.size(); ++i ) 00195 nodes.insert( connect[i] ); 00196 } 00197 } 00198 00199 if( nodes.empty() ) 00200 { 00201 MB_SET_ERR( MB_ENTITY_NOT_FOUND, "Nothing to write" ); 00202 } 00203 00204 return MB_SUCCESS; 00205 } 00206 00207 ErrorCode WriteVtk::write_header( std::ostream& stream ) 00208 { 00209 stream << "# vtk DataFile Version 3.0" << std::endl; 00210 stream << MOAB_VERSION_STRING << std::endl; 00211 stream << "ASCII" << std::endl; 00212 stream << "DATASET UNSTRUCTURED_GRID" << std::endl; 00213 return MB_SUCCESS; 00214 } 00215 00216 ErrorCode WriteVtk::write_nodes( std::ostream& stream, const Range& nodes ) 00217 { 00218 ErrorCode rval; 00219 00220 stream << "POINTS " << nodes.size() << " double" << std::endl; 00221 00222 double coords[3]; 00223 for( Range::const_iterator i = nodes.begin(); i != nodes.end(); ++i ) 00224 { 00225 coords[1] = coords[2] = 0.0; 00226 rval = mbImpl->get_coords( &( *i ), 1, coords ); 00227 if( MB_SUCCESS != rval ) return rval; 00228 stream << coords[0] << ' ' << coords[1] << ' ' << coords[2] << std::endl; 00229 } 00230 00231 return MB_SUCCESS; 00232 } 00233 00234 ErrorCode WriteVtk::write_elems( std::ostream& stream, const Range& nodes, const Range& elems ) 00235 { 00236 ErrorCode rval; 00237 00238 Range connectivity; // because we now support polyhedra, it could contain faces 00239 rval = mbImpl->get_connectivity( elems, connectivity );MB_CHK_ERR( rval ); 00240 00241 Range nodes_from_connectivity = connectivity.subset_by_type( MBVERTEX ); 00242 Range faces_from_connectivity = 00243 subtract( connectivity, nodes_from_connectivity ); // these could be faces of polyhedra 00244 00245 Range connected_nodes; 00246 rval = mbImpl->get_connectivity( faces_from_connectivity, connected_nodes );MB_CHK_ERR( rval ); 00247 connected_nodes.merge( nodes_from_connectivity ); 00248 00249 Range free_nodes = subtract( nodes, connected_nodes ); 00250 00251 // Get and write counts 00252 unsigned long num_elems, num_uses; 00253 num_elems = num_uses = elems.size(); 00254 00255 std::map< EntityHandle, int > sizeFieldsPolyhedra; 00256 00257 for( Range::const_iterator i = elems.begin(); i != elems.end(); ++i ) 00258 { 00259 EntityType type = mbImpl->type_from_handle( *i ); 00260 if( !VtkUtil::get_vtk_type( type, CN::VerticesPerEntity( type ) ) ) continue; 00261 00262 EntityHandle elem = *i; 00263 const EntityHandle* connect = NULL; 00264 int conn_len = 0; 00265 // Dummy storage vector for structured mesh "get_connectivity" function 00266 std::vector< EntityHandle > storage; 00267 rval = mbImpl->get_connectivity( elem, connect, conn_len, false, &storage );MB_CHK_ERR( rval ); 00268 00269 num_uses += conn_len; 00270 // if polyhedra, we will count the number of nodes in each face too 00271 if( TYPE_FROM_HANDLE( elem ) == MBPOLYHEDRON ) 00272 { 00273 int numFields = 1; // there will be one for number of faces; forgot about this one 00274 for( int j = 0; j < conn_len; j++ ) 00275 { 00276 const EntityHandle* conn = NULL; 00277 int num_nd = 0; 00278 rval = mbImpl->get_connectivity( connect[j], conn, num_nd );MB_CHK_ERR( rval ); 00279 numFields += num_nd + 1; 00280 } 00281 sizeFieldsPolyhedra[elem] = numFields; // will be used later, at writing 00282 num_uses += ( numFields - conn_len ); 00283 } 00284 } 00285 freeNodes = (int)free_nodes.size(); 00286 if( !createOneNodeCells ) freeNodes = 0; // do not create one node cells 00287 stream << "CELLS " << num_elems + freeNodes << ' ' << num_uses + 2 * freeNodes << std::endl; 00288 00289 // Write element connectivity 00290 std::vector< int > conn_data; 00291 std::vector< unsigned > vtk_types( elems.size() + freeNodes ); 00292 std::vector< unsigned >::iterator t = vtk_types.begin(); 00293 for( Range::const_iterator i = elems.begin(); i != elems.end(); ++i ) 00294 { 00295 // Get type information for element 00296 EntityHandle elem = *i; 00297 EntityType type = TYPE_FROM_HANDLE( elem ); 00298 00299 // Get element connectivity 00300 const EntityHandle* connect = NULL; 00301 int conn_len = 0; 00302 // Dummy storage vector for structured mesh "get_connectivity" function 00303 std::vector< EntityHandle > storage; 00304 rval = mbImpl->get_connectivity( elem, connect, conn_len, false, &storage );MB_CHK_ERR( rval ); 00305 00306 // Get VTK type 00307 const VtkElemType* vtk_type = VtkUtil::get_vtk_type( type, conn_len ); 00308 if( !vtk_type ) 00309 { 00310 // Try connectivity with 1 fewer node 00311 vtk_type = VtkUtil::get_vtk_type( type, conn_len - 1 ); 00312 if( vtk_type ) 00313 conn_len--; 00314 else 00315 { 00316 MB_SET_ERR( MB_FAILURE, "Vtk file format does not support elements of type " 00317 << CN::EntityTypeName( type ) << " (" << (int)type << ") with " << conn_len 00318 << " nodes" ); 00319 } 00320 } 00321 00322 // Save VTK type index for later 00323 *t = vtk_type->vtk_type; 00324 ++t; 00325 00326 if( type != MBPOLYHEDRON ) 00327 { 00328 // Get IDs from vertex handles 00329 assert( conn_len > 0 ); 00330 conn_data.resize( conn_len ); 00331 for( int j = 0; j < conn_len; ++j ) 00332 conn_data[j] = nodes.index( connect[j] ); 00333 00334 // Write connectivity list 00335 stream << conn_len; 00336 if( vtk_type->node_order ) 00337 for( int k = 0; k < conn_len; ++k ) 00338 stream << ' ' << conn_data[vtk_type->node_order[k]]; 00339 else 00340 for( int k = 0; k < conn_len; ++k ) 00341 stream << ' ' << conn_data[k]; 00342 stream << std::endl; 00343 } 00344 else 00345 { 00346 // POLYHEDRON needs a special case, loop over faces to get nodes 00347 stream << sizeFieldsPolyhedra[elem] << " " << conn_len; 00348 for( int k = 0; k < conn_len; k++ ) 00349 { 00350 EntityHandle face = connect[k]; 00351 const EntityHandle* conn = NULL; 00352 int num_nodes = 0; 00353 rval = mbImpl->get_connectivity( face, conn, num_nodes );MB_CHK_ERR( rval ); 00354 // num_uses += num_nd + 1; // 1 for number of vertices in face 00355 conn_data.resize( num_nodes ); 00356 for( int j = 0; j < num_nodes; ++j ) 00357 conn_data[j] = nodes.index( conn[j] ); 00358 00359 stream << ' ' << num_nodes; 00360 00361 for( int j = 0; j < num_nodes; ++j ) 00362 stream << ' ' << conn_data[j]; 00363 } 00364 stream << std::endl; 00365 } 00366 } 00367 00368 if( createOneNodeCells ) 00369 for( Range::const_iterator v = free_nodes.begin(); v != free_nodes.end(); ++v, ++t ) 00370 { 00371 EntityHandle node = *v; 00372 stream << "1 " << nodes.index( node ) << std::endl; 00373 *t = 1; 00374 } 00375 00376 // Write element types 00377 stream << "CELL_TYPES " << vtk_types.size() << std::endl; 00378 for( std::vector< unsigned >::const_iterator i = vtk_types.begin(); i != vtk_types.end(); ++i ) 00379 stream << *i << std::endl; 00380 00381 return MB_SUCCESS; 00382 } 00383 00384 ErrorCode WriteVtk::write_tags( std::ostream& stream, 00385 bool nodes, 00386 const Range& entities, 00387 const Tag* tag_list, 00388 int num_tags ) 00389 { 00390 ErrorCode rval; 00391 00392 // The #$%@#$% MOAB interface does not have a function to retrieve 00393 // all entities with a tag, only all entities with a specified type 00394 // and tag. Define types to loop over depending on the if vertex 00395 // or element tag data is being written. It seems horribly inefficient 00396 // to have the implementation subdivide the results by type and have 00397 // to call that function once for each type just to recombine the results. 00398 // Unfortunately, there doesn't seem to be any other way. 00399 EntityType low_type, high_type; 00400 if( nodes ) 00401 { 00402 low_type = MBVERTEX; 00403 high_type = MBEDGE; 00404 } 00405 else 00406 { 00407 low_type = MBEDGE; 00408 high_type = MBENTITYSET; 00409 } 00410 00411 // Get all defined tags 00412 std::vector< Tag > tags; 00413 std::vector< Tag >::iterator i; 00414 rval = writeTool->get_tag_list( tags, tag_list, num_tags, false ); 00415 if( MB_SUCCESS != rval ) return rval; 00416 00417 // For each tag... 00418 bool entities_have_tags = false; 00419 for( i = tags.begin(); i != tags.end(); ++i ) 00420 { 00421 // Skip tags holding entity handles -- no way to save them 00422 DataType dtype; 00423 rval = mbImpl->tag_get_data_type( *i, dtype ); 00424 if( MB_SUCCESS != rval ) return rval; 00425 if( dtype == MB_TYPE_HANDLE ) continue; 00426 00427 // If in strict mode, don't write tags that do not fit in any 00428 // attribute type (SCALAR : 1 to 4 values, VECTOR : 3 values, TENSOR : 9 values) 00429 if( mStrict ) 00430 { 00431 int count; 00432 rval = mbImpl->tag_get_length( *i, count ); 00433 if( MB_SUCCESS != rval ) return rval; 00434 if( count < 1 || ( count > 4 && count != 9 ) ) continue; 00435 } 00436 00437 // Get subset of input entities that have the tag set 00438 Range tagged; 00439 for( EntityType type = low_type; type < high_type; ++type ) 00440 { 00441 Range tmp_tagged; 00442 rval = mbImpl->get_entities_by_type_and_tag( 0, type, &( *i ), 0, 1, tmp_tagged ); 00443 if( MB_SUCCESS != rval ) return rval; 00444 tmp_tagged = intersect( tmp_tagged, entities ); 00445 tagged.merge( tmp_tagged ); 00446 } 00447 00448 // If any entities were tagged 00449 if( !tagged.empty() ) 00450 { 00451 // If this is the first tag being written for the 00452 // entity type, write the label marking the beginning 00453 // of the tag data. 00454 if( !entities_have_tags ) 00455 { 00456 entities_have_tags = true; 00457 if( nodes ) 00458 stream << "POINT_DATA " << entities.size() << std::endl; 00459 else 00460 stream << "CELL_DATA " << entities.size() + freeNodes << std::endl; 00461 } 00462 00463 // Write the tag 00464 rval = write_tag( stream, *i, entities, tagged ); 00465 if( MB_SUCCESS != rval ) return rval; 00466 } 00467 } 00468 00469 return MB_SUCCESS; 00470 } 00471 00472 template < typename T > 00473 void WriteVtk::write_data( std::ostream& stream, const std::vector< T >& data, unsigned vals_per_tag ) 00474 { 00475 typename std::vector< T >::const_iterator d = data.begin(); 00476 const unsigned long n = data.size() / vals_per_tag; 00477 00478 for( unsigned long i = 0; i < n; ++i ) 00479 { 00480 for( unsigned j = 0; j < vals_per_tag; ++j, ++d ) 00481 { 00482 if( sizeof( T ) == 1 ) 00483 stream << (unsigned int)*d << ' '; 00484 else 00485 stream << *d << ' '; 00486 } 00487 stream << std::endl; 00488 } 00489 } 00490 00491 // template <> 00492 // void WriteVtk::write_data<unsigned char>(std::ostream& stream, 00493 // const std::vector<unsigned char>& data, 00494 // unsigned vals_per_tag) 00495 //{ 00496 // std::vector<unsigned char>::const_iterator d = data.begin(); 00497 // const unsigned long n = data.size() / vals_per_tag; 00498 // 00499 // for (unsigned long i = 0; i < n; ++i) { 00500 // for (unsigned j = 0; j < vals_per_tag; ++j, ++d) 00501 // stream << (unsigned int)*d << ' '; 00502 // stream << std::endl; 00503 // } 00504 //} 00505 00506 template < typename T > 00507 ErrorCode WriteVtk::write_tag( std::ostream& stream, Tag tag, const Range& entities, const Range& tagged, const int ) 00508 { 00509 ErrorCode rval; 00510 int addFreeNodes = 0; 00511 if( TYPE_FROM_HANDLE( entities[0] ) > MBVERTEX ) addFreeNodes = freeNodes; 00512 // we created freeNodes 1-node cells, so we have to augment cell data too 00513 // we know that the 1 node cells are added at the end, after all other cells; 00514 // so the default values will be set to those extra , artificial cells 00515 const unsigned long n = entities.size() + addFreeNodes; 00516 00517 // Get tag properties 00518 00519 std::string name; 00520 int vals_per_tag; 00521 if( MB_SUCCESS != mbImpl->tag_get_name( tag, name ) || MB_SUCCESS != mbImpl->tag_get_length( tag, vals_per_tag ) ) 00522 return MB_FAILURE; 00523 00524 // Get a tag value for each entity. Do this by initializing the 00525 // "data" vector with zero, and then filling in the values for 00526 // the entities that actually have the tag set. 00527 std::vector< T > data; 00528 data.resize( n * vals_per_tag, 0 ); 00529 // If there is a default value for the tag, set the actual default value 00530 std::vector< T > def_value( vals_per_tag ); 00531 rval = mbImpl->tag_get_default_value( tag, &( def_value[0] ) ); 00532 if( MB_SUCCESS == rval ) SysUtil::setmem( &( data[0] ), &( def_value[0] ), vals_per_tag * sizeof( T ), n ); 00533 00534 Range::const_iterator t = tagged.begin(); 00535 typename std::vector< T >::iterator d = data.begin(); 00536 for( Range::const_iterator i = entities.begin(); i != entities.end() && t != tagged.end(); ++i, d += vals_per_tag ) 00537 { 00538 if( *i == *t ) 00539 { 00540 ++t; 00541 rval = mbImpl->tag_get_data( tag, &( *i ), 1, &( *d ) ); 00542 if( MB_SUCCESS != rval ) return rval; 00543 } 00544 } 00545 00546 // Write the tag values, one entity per line. 00547 write_data( stream, data, vals_per_tag ); 00548 00549 return MB_SUCCESS; 00550 } 00551 00552 ErrorCode WriteVtk::write_bit_tag( std::ostream& stream, Tag tag, const Range& entities, const Range& tagged ) 00553 { 00554 ErrorCode rval; 00555 const unsigned long n = entities.size(); 00556 00557 // Get tag properties 00558 00559 std::string name; 00560 int vals_per_tag; 00561 if( MB_SUCCESS != mbImpl->tag_get_name( tag, name ) || MB_SUCCESS != mbImpl->tag_get_length( tag, vals_per_tag ) ) 00562 return MB_FAILURE; 00563 00564 if( vals_per_tag > 8 ) 00565 { 00566 MB_SET_ERR( MB_FAILURE, "Invalid tag size for bit tag \"" << name << "\"" ); 00567 } 00568 00569 // Get a tag value for each entity. 00570 // Get bits for each entity and unpack into 00571 // one integer in the 'data' array for each bit. 00572 // Initialize 'data' to zero because we will skip 00573 // those entities for which the tag is not set. 00574 std::vector< unsigned short > data; 00575 data.resize( n * vals_per_tag, 0 ); 00576 Range::const_iterator t = tagged.begin(); 00577 std::vector< unsigned short >::iterator d = data.begin(); 00578 for( Range::const_iterator i = entities.begin(); i != entities.end() && t != tagged.end(); ++i ) 00579 { 00580 if( *i == *t ) 00581 { 00582 ++t; 00583 unsigned char value; 00584 rval = mbImpl->tag_get_data( tag, &( *i ), 1, &value ); 00585 for( int j = 0; j < vals_per_tag; ++j, ++d ) 00586 *d = (unsigned short)( value & ( 1 << j ) ? 1 : 0 ); 00587 if( MB_SUCCESS != rval ) return rval; 00588 } 00589 else 00590 { 00591 // If tag is not set for entity, skip values in array 00592 d += vals_per_tag; 00593 } 00594 } 00595 00596 // Write the tag values, one entity per line. 00597 write_data( stream, data, vals_per_tag ); 00598 00599 return MB_SUCCESS; 00600 } 00601 00602 ErrorCode WriteVtk::write_tag( std::ostream& s, Tag tag, const Range& entities, const Range& tagged ) 00603 { 00604 // Get tag properties 00605 std::string name; 00606 DataType type; 00607 int vals_per_tag; 00608 if( MB_SUCCESS != mbImpl->tag_get_name( tag, name ) || MB_SUCCESS != mbImpl->tag_get_length( tag, vals_per_tag ) || 00609 MB_SUCCESS != mbImpl->tag_get_data_type( tag, type ) ) 00610 return MB_FAILURE; 00611 00612 // Skip tags of type ENTITY_HANDLE 00613 if( MB_TYPE_HANDLE == type ) return MB_FAILURE; 00614 00615 // Now that we're past the point where the name would be used in 00616 // an error message, remove any spaces to conform to VTK file. 00617 for( std::string::iterator i = name.begin(); i != name.end(); ++i ) 00618 { 00619 if( isspace( *i ) || iscntrl( *i ) ) *i = '_'; 00620 } 00621 00622 // Write the tag description 00623 if( 3 == vals_per_tag && MB_TYPE_DOUBLE == type ) 00624 s << "VECTORS " << name << ' ' << VtkUtil::vtkTypeNames[type] << std::endl; 00625 else if( 9 == vals_per_tag ) 00626 s << "TENSORS " << name << ' ' << VtkUtil::vtkTypeNames[type] << std::endl; 00627 else 00628 s << "SCALARS " << name << ' ' << VtkUtil::vtkTypeNames[type] << ' ' << vals_per_tag << std::endl 00629 << "LOOKUP_TABLE default" << std::endl; 00630 00631 // Write the tag data 00632 switch( type ) 00633 { 00634 case MB_TYPE_OPAQUE: 00635 return write_tag< unsigned char >( s, tag, entities, tagged, 0 ); 00636 case MB_TYPE_INTEGER: 00637 return write_tag< int >( s, tag, entities, tagged, 0 ); 00638 case MB_TYPE_DOUBLE: 00639 return write_tag< double >( s, tag, entities, tagged, 0 ); 00640 case MB_TYPE_BIT: 00641 return write_bit_tag( s, tag, entities, tagged ); 00642 default: 00643 return MB_FAILURE; 00644 } 00645 } 00646 00647 } // namespace moab