Mesh Oriented datABase  (version 5.5.0)
An array-based unstructured mesh library
WriteVtk.cpp
Go to the documentation of this file.
1 /**
2  * MOAB, a Mesh-Oriented datABase, is a software component for creating,
3  * storing and accessing finite element mesh data.
4  *
5  * Copyright 2004 Sandia Corporation. Under the terms of Contract
6  * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government
7  * retains certain rights in this software.
8  *
9  * This library is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public
11  * License as published by the Free Software Foundation; either
12  * version 2.1 of the License, or (at your option) any later version.
13  *
14  */
15 
16 #ifdef WIN32
17 #ifdef _DEBUG
18 // turn off warnings that say they debugging identifier has been truncated
19 // this warning comes up when using some STL containers
20 #pragma warning( disable : 4786 )
21 #endif
22 #endif
23 
24 #include "WriteVtk.hpp"
25 #include "moab/VtkUtil.hpp"
26 #include "SysUtil.hpp"
27 
28 #include <fstream>
29 #include <iostream>
30 #include <cstdio>
31 #include <cassert>
32 #include <vector>
33 #include <set>
34 #include <map>
35 #include <iterator>
36 
37 #include "moab/Interface.hpp"
38 #include "moab/Range.hpp"
39 #include "moab/CN.hpp"
40 #include "MBTagConventions.hpp"
41 #include "moab/WriteUtilIface.hpp"
42 #include "Internals.hpp"
43 #include "moab/FileOptions.hpp"
44 
45 #define INS_ID( stringvar, prefix, id ) sprintf( stringvar, prefix, id )
46 
47 namespace moab
48 {
49 
50 const int DEFAULT_PRECISION = 10;
51 const bool DEFAULT_STRICT = true;
52 
54 {
55  return new WriteVtk( iface );
56 }
57 
59  : mbImpl( impl ), writeTool( 0 ), mStrict( DEFAULT_STRICT ), freeNodes( 0 ), createOneNodeCells( false )
60 {
61  assert( impl != NULL );
62  impl->query_interface( writeTool );
63 }
64 
66 {
68 }
69 
70 ErrorCode WriteVtk::write_file( const char* file_name,
71  const bool overwrite,
72  const FileOptions& opts,
73  const EntityHandle* output_list,
74  const int num_sets,
75  const std::vector< std::string >& /* qa_list */,
76  const Tag* tag_list,
77  int num_tags,
78  int /* export_dimension */ )
79 {
80  ErrorCode rval;
81 
82  // Get precision for node coordinates
83  int precision;
84  if( MB_SUCCESS != opts.get_int_option( "PRECISION", precision ) ) precision = DEFAULT_PRECISION;
85 
86  if( MB_SUCCESS == opts.get_null_option( "STRICT" ) )
87  mStrict = true;
88  else if( MB_SUCCESS == opts.get_null_option( "RELAXED" ) )
89  mStrict = false;
90  else
92 
93  if( MB_SUCCESS == opts.get_null_option( "CREATE_ONE_NODE_CELLS" ) ) createOneNodeCells = true;
94 
95  // Get entities to write
96  Range nodes, elems;
97  rval = gather_mesh( output_list, num_sets, nodes, elems );
98  if( MB_SUCCESS != rval ) return rval;
99 
100  // Honor overwrite flag
101  if( !overwrite )
102  {
103  rval = writeTool->check_doesnt_exist( file_name );
104  if( MB_SUCCESS != rval ) return rval;
105  }
106 
107  // Create file
108  std::ofstream file( file_name );
109  if( !file )
110  {
111  MB_SET_ERR( MB_FILE_WRITE_ERROR, "Could not open file: " << file_name );
112  }
113  file.precision( precision );
114 
115  // Write file
116  if( ( rval = write_header( file ) ) != MB_SUCCESS || ( rval = write_nodes( file, nodes ) ) != MB_SUCCESS ||
117  ( rval = write_elems( file, nodes, elems ) ) != MB_SUCCESS ||
118  ( rval = write_tags( file, true, nodes, tag_list, num_tags ) ) != MB_SUCCESS ||
119  ( rval = write_tags( file, false, elems, tag_list, num_tags ) ) != MB_SUCCESS )
120  {
121  file.close();
122  remove( file_name );
123  return rval;
124  }
125 
126  return MB_SUCCESS;
127 }
128 
129 ErrorCode WriteVtk::gather_mesh( const EntityHandle* set_list, int num_sets, Range& nodes, Range& elems )
130 {
131  ErrorCode rval;
132  int e;
133 
134  if( !set_list || !num_sets )
135  {
136  Range a;
137  rval = mbImpl->get_entities_by_handle( 0, a );
138  if( MB_SUCCESS != rval ) return rval;
139 
140  Range::const_iterator node_i, elem_i, set_i;
141  node_i = a.lower_bound( a.begin(), a.end(), CREATE_HANDLE( MBVERTEX, 0, e ) );
142  elem_i = a.lower_bound( node_i, a.end(), CREATE_HANDLE( MBEDGE, 0, e ) );
143  set_i = a.lower_bound( elem_i, a.end(), CREATE_HANDLE( MBENTITYSET, 0, e ) );
144  nodes.merge( node_i, elem_i );
145  elems.merge( elem_i, set_i );
146 
147  // Filter out unsupported element types
148  EntityType et = MBEDGE;
149  for( et++; et < MBENTITYSET; et++ )
150  {
151  if( VtkUtil::get_vtk_type( et, CN::VerticesPerEntity( et ) ) ) continue;
152  Range::iterator eit = elems.lower_bound( elems.begin(), elems.end(), CREATE_HANDLE( et, 0, e ) ),
153  ep1it = elems.lower_bound( elems.begin(), elems.end(), CREATE_HANDLE( et + 1, 0, e ) );
154  elems.erase( eit, ep1it );
155  }
156  }
157  else
158  {
159  std::set< EntityHandle > visited;
160  std::vector< EntityHandle > sets;
161  sets.reserve( num_sets );
162  std::copy( set_list, set_list + num_sets, std::back_inserter( sets ) );
163  while( !sets.empty() )
164  {
165  // Get next set
166  EntityHandle set = sets.back();
167  sets.pop_back();
168  // Skip sets we've already done
169  if( !visited.insert( set ).second ) continue;
170 
171  Range a;
172  rval = mbImpl->get_entities_by_handle( set, a );
173  if( MB_SUCCESS != rval ) return rval;
174 
175  Range::const_iterator node_i, elem_i, set_i;
176  node_i = a.lower_bound( a.begin(), a.end(), CREATE_HANDLE( MBVERTEX, 0, e ) );
177  elem_i = a.lower_bound( node_i, a.end(), CREATE_HANDLE( MBEDGE, 0, e ) );
178  set_i = a.lower_bound( elem_i, a.end(), CREATE_HANDLE( MBENTITYSET, 0, e ) );
179  nodes.merge( node_i, elem_i );
180  elems.merge( elem_i, set_i );
181  std::copy( set_i, a.end(), std::back_inserter( sets ) );
182 
183  a.clear();
184  rval = mbImpl->get_child_meshsets( set, a );
185  std::copy( a.begin(), a.end(), std::back_inserter( sets ) );
186  }
187 
188  for( Range::const_iterator ei = elems.begin(); ei != elems.end(); ++ei )
189  {
190  std::vector< EntityHandle > connect;
191  rval = mbImpl->get_connectivity( &( *ei ), 1, connect );
192  if( MB_SUCCESS != rval ) return rval;
193 
194  for( unsigned int i = 0; i < connect.size(); ++i )
195  nodes.insert( connect[i] );
196  }
197  }
198 
199  if( nodes.empty() )
200  {
201  MB_SET_ERR( MB_ENTITY_NOT_FOUND, "Nothing to write" );
202  }
203 
204  return MB_SUCCESS;
205 }
206 
207 ErrorCode WriteVtk::write_header( std::ostream& stream )
208 {
209  stream << "# vtk DataFile Version 3.0" << std::endl;
210  stream << MOAB_VERSION_STRING << std::endl;
211  stream << "ASCII" << std::endl;
212  stream << "DATASET UNSTRUCTURED_GRID" << std::endl;
213  return MB_SUCCESS;
214 }
215 
216 ErrorCode WriteVtk::write_nodes( std::ostream& stream, const Range& nodes )
217 {
218  ErrorCode rval;
219 
220  stream << "POINTS " << nodes.size() << " double" << std::endl;
221 
222  double coords[3];
223  for( Range::const_iterator i = nodes.begin(); i != nodes.end(); ++i )
224  {
225  coords[1] = coords[2] = 0.0;
226  rval = mbImpl->get_coords( &( *i ), 1, coords );
227  if( MB_SUCCESS != rval ) return rval;
228  stream << coords[0] << ' ' << coords[1] << ' ' << coords[2] << std::endl;
229  }
230 
231  return MB_SUCCESS;
232 }
233 
234 ErrorCode WriteVtk::write_elems( std::ostream& stream, const Range& nodes, const Range& elems )
235 {
236  ErrorCode rval;
237 
238  Range connectivity; // because we now support polyhedra, it could contain faces
239  rval = mbImpl->get_connectivity( elems, connectivity );MB_CHK_ERR( rval );
240 
241  Range nodes_from_connectivity = connectivity.subset_by_type( MBVERTEX );
242  Range faces_from_connectivity =
243  subtract( connectivity, nodes_from_connectivity ); // these could be faces of polyhedra
244 
245  Range connected_nodes;
246  rval = mbImpl->get_connectivity( faces_from_connectivity, connected_nodes );MB_CHK_ERR( rval );
247  connected_nodes.merge( nodes_from_connectivity );
248 
249  Range free_nodes = subtract( nodes, connected_nodes );
250 
251  // Get and write counts
252  unsigned long num_elems, num_uses;
253  num_elems = num_uses = elems.size();
254 
255  std::map< EntityHandle, int > sizeFieldsPolyhedra;
256 
257  for( Range::const_iterator i = elems.begin(); i != elems.end(); ++i )
258  {
259  EntityType type = mbImpl->type_from_handle( *i );
260  if( !VtkUtil::get_vtk_type( type, CN::VerticesPerEntity( type ) ) ) continue;
261 
262  EntityHandle elem = *i;
263  const EntityHandle* connect = NULL;
264  int conn_len = 0;
265  // Dummy storage vector for structured mesh "get_connectivity" function
266  std::vector< EntityHandle > storage;
267  rval = mbImpl->get_connectivity( elem, connect, conn_len, false, &storage );MB_CHK_ERR( rval );
268 
269  num_uses += conn_len;
270  // if polyhedra, we will count the number of nodes in each face too
271  if( TYPE_FROM_HANDLE( elem ) == MBPOLYHEDRON )
272  {
273  int numFields = 1; // there will be one for number of faces; forgot about this one
274  for( int j = 0; j < conn_len; j++ )
275  {
276  const EntityHandle* conn = NULL;
277  int num_nd = 0;
278  rval = mbImpl->get_connectivity( connect[j], conn, num_nd );MB_CHK_ERR( rval );
279  numFields += num_nd + 1;
280  }
281  sizeFieldsPolyhedra[elem] = numFields; // will be used later, at writing
282  num_uses += ( numFields - conn_len );
283  }
284  }
285  freeNodes = (int)free_nodes.size();
286  if( !createOneNodeCells ) freeNodes = 0; // do not create one node cells
287  stream << "CELLS " << num_elems + freeNodes << ' ' << num_uses + 2 * freeNodes << std::endl;
288 
289  // Write element connectivity
290  std::vector< int > conn_data;
291  std::vector< unsigned > vtk_types( elems.size() + freeNodes );
292  std::vector< unsigned >::iterator t = vtk_types.begin();
293  for( Range::const_iterator i = elems.begin(); i != elems.end(); ++i )
294  {
295  // Get type information for element
296  EntityHandle elem = *i;
297  EntityType type = TYPE_FROM_HANDLE( elem );
298 
299  // Get element connectivity
300  const EntityHandle* connect = NULL;
301  int conn_len = 0;
302  // Dummy storage vector for structured mesh "get_connectivity" function
303  std::vector< EntityHandle > storage;
304  rval = mbImpl->get_connectivity( elem, connect, conn_len, false, &storage );MB_CHK_ERR( rval );
305 
306  // Get VTK type
307  const VtkElemType* vtk_type = VtkUtil::get_vtk_type( type, conn_len );
308  if( !vtk_type )
309  {
310  // Try connectivity with 1 fewer node
311  vtk_type = VtkUtil::get_vtk_type( type, conn_len - 1 );
312  if( vtk_type )
313  conn_len--;
314  else
315  {
316  MB_SET_ERR( MB_FAILURE, "Vtk file format does not support elements of type "
317  << CN::EntityTypeName( type ) << " (" << (int)type << ") with " << conn_len
318  << " nodes" );
319  }
320  }
321 
322  // Save VTK type index for later
323  *t = vtk_type->vtk_type;
324  ++t;
325 
326  if( type != MBPOLYHEDRON )
327  {
328  // Get IDs from vertex handles
329  assert( conn_len > 0 );
330  conn_data.resize( conn_len );
331  for( int j = 0; j < conn_len; ++j )
332  conn_data[j] = nodes.index( connect[j] );
333 
334  // Write connectivity list
335  stream << conn_len;
336  if( vtk_type->node_order )
337  for( int k = 0; k < conn_len; ++k )
338  stream << ' ' << conn_data[vtk_type->node_order[k]];
339  else
340  for( int k = 0; k < conn_len; ++k )
341  stream << ' ' << conn_data[k];
342  stream << std::endl;
343  }
344  else
345  {
346  // POLYHEDRON needs a special case, loop over faces to get nodes
347  stream << sizeFieldsPolyhedra[elem] << " " << conn_len;
348  for( int k = 0; k < conn_len; k++ )
349  {
350  EntityHandle face = connect[k];
351  const EntityHandle* conn = NULL;
352  int num_nodes = 0;
353  rval = mbImpl->get_connectivity( face, conn, num_nodes );MB_CHK_ERR( rval );
354  // num_uses += num_nd + 1; // 1 for number of vertices in face
355  conn_data.resize( num_nodes );
356  for( int j = 0; j < num_nodes; ++j )
357  conn_data[j] = nodes.index( conn[j] );
358 
359  stream << ' ' << num_nodes;
360 
361  for( int j = 0; j < num_nodes; ++j )
362  stream << ' ' << conn_data[j];
363  }
364  stream << std::endl;
365  }
366  }
367 
368  if( createOneNodeCells )
369  for( Range::const_iterator v = free_nodes.begin(); v != free_nodes.end(); ++v, ++t )
370  {
371  EntityHandle node = *v;
372  stream << "1 " << nodes.index( node ) << std::endl;
373  *t = 1;
374  }
375 
376  // Write element types
377  stream << "CELL_TYPES " << vtk_types.size() << std::endl;
378  for( std::vector< unsigned >::const_iterator i = vtk_types.begin(); i != vtk_types.end(); ++i )
379  stream << *i << std::endl;
380 
381  return MB_SUCCESS;
382 }
383 
384 ErrorCode WriteVtk::write_tags( std::ostream& stream,
385  bool nodes,
386  const Range& entities,
387  const Tag* tag_list,
388  int num_tags )
389 {
390  ErrorCode rval;
391 
392  // The #$%@#$% MOAB interface does not have a function to retrieve
393  // all entities with a tag, only all entities with a specified type
394  // and tag. Define types to loop over depending on the if vertex
395  // or element tag data is being written. It seems horribly inefficient
396  // to have the implementation subdivide the results by type and have
397  // to call that function once for each type just to recombine the results.
398  // Unfortunately, there doesn't seem to be any other way.
399  EntityType low_type, high_type;
400  if( nodes )
401  {
402  low_type = MBVERTEX;
403  high_type = MBEDGE;
404  }
405  else
406  {
407  low_type = MBEDGE;
408  high_type = MBENTITYSET;
409  }
410 
411  // Get all defined tags
412  std::vector< Tag > tags;
413  std::vector< Tag >::iterator i;
414  rval = writeTool->get_tag_list( tags, tag_list, num_tags, false );
415  if( MB_SUCCESS != rval ) return rval;
416 
417  // For each tag...
418  bool entities_have_tags = false;
419  for( i = tags.begin(); i != tags.end(); ++i )
420  {
421  // Skip tags holding entity handles -- no way to save them
422  DataType dtype;
423  rval = mbImpl->tag_get_data_type( *i, dtype );
424  if( MB_SUCCESS != rval ) return rval;
425  if( dtype == MB_TYPE_HANDLE ) continue;
426 
427  // If in strict mode, don't write tags that do not fit in any
428  // attribute type (SCALAR : 1 to 4 values, VECTOR : 3 values, TENSOR : 9 values)
429  if( mStrict )
430  {
431  int count;
432  rval = mbImpl->tag_get_length( *i, count );
433  if( MB_SUCCESS != rval ) return rval;
434  if( count < 1 || ( count > 4 && count != 9 ) ) continue;
435  }
436 
437  // Get subset of input entities that have the tag set
438  Range tagged;
439  for( EntityType type = low_type; type < high_type; ++type )
440  {
441  Range tmp_tagged;
442  rval = mbImpl->get_entities_by_type_and_tag( 0, type, &( *i ), 0, 1, tmp_tagged );
443  if( MB_SUCCESS != rval ) return rval;
444  tmp_tagged = intersect( tmp_tagged, entities );
445  tagged.merge( tmp_tagged );
446  }
447 
448  // If any entities were tagged
449  if( !tagged.empty() )
450  {
451  // If this is the first tag being written for the
452  // entity type, write the label marking the beginning
453  // of the tag data.
454  if( !entities_have_tags )
455  {
456  entities_have_tags = true;
457  if( nodes )
458  stream << "POINT_DATA " << entities.size() << std::endl;
459  else
460  stream << "CELL_DATA " << entities.size() + freeNodes << std::endl;
461  }
462 
463  // Write the tag
464  rval = write_tag( stream, *i, entities, tagged );
465  if( MB_SUCCESS != rval ) return rval;
466  }
467  }
468 
469  return MB_SUCCESS;
470 }
471 
472 template < typename T >
473 void WriteVtk::write_data( std::ostream& stream, const std::vector< T >& data, unsigned vals_per_tag )
474 {
475  typename std::vector< T >::const_iterator d = data.begin();
476  const unsigned long n = data.size() / vals_per_tag;
477 
478  for( unsigned long i = 0; i < n; ++i )
479  {
480  for( unsigned j = 0; j < vals_per_tag; ++j, ++d )
481  {
482  if( sizeof( T ) == 1 )
483  stream << (unsigned int)*d << ' ';
484  else
485  stream << *d << ' ';
486  }
487  stream << std::endl;
488  }
489 }
490 
491 // template <>
492 // void WriteVtk::write_data<unsigned char>(std::ostream& stream,
493 // const std::vector<unsigned char>& data,
494 // unsigned vals_per_tag)
495 //{
496 // std::vector<unsigned char>::const_iterator d = data.begin();
497 // const unsigned long n = data.size() / vals_per_tag;
498 //
499 // for (unsigned long i = 0; i < n; ++i) {
500 // for (unsigned j = 0; j < vals_per_tag; ++j, ++d)
501 // stream << (unsigned int)*d << ' ';
502 // stream << std::endl;
503 // }
504 //}
505 
506 template < typename T >
507 ErrorCode WriteVtk::write_tag( std::ostream& stream, Tag tag, const Range& entities, const Range& tagged, const int )
508 {
509  ErrorCode rval;
510  int addFreeNodes = 0;
511  if( TYPE_FROM_HANDLE( entities[0] ) > MBVERTEX ) addFreeNodes = freeNodes;
512  // we created freeNodes 1-node cells, so we have to augment cell data too
513  // we know that the 1 node cells are added at the end, after all other cells;
514  // so the default values will be set to those extra , artificial cells
515  const unsigned long n = entities.size() + addFreeNodes;
516 
517  // Get tag properties
518 
519  std::string name;
520  int vals_per_tag;
521  if( MB_SUCCESS != mbImpl->tag_get_name( tag, name ) || MB_SUCCESS != mbImpl->tag_get_length( tag, vals_per_tag ) )
522  return MB_FAILURE;
523 
524  // Get a tag value for each entity. Do this by initializing the
525  // "data" vector with zero, and then filling in the values for
526  // the entities that actually have the tag set.
527  std::vector< T > data;
528  data.resize( n * vals_per_tag, 0 );
529  // If there is a default value for the tag, set the actual default value
530  std::vector< T > def_value( vals_per_tag );
531  rval = mbImpl->tag_get_default_value( tag, &( def_value[0] ) );
532  if( MB_SUCCESS == rval ) SysUtil::setmem( &( data[0] ), &( def_value[0] ), vals_per_tag * sizeof( T ), n );
533 
534  Range::const_iterator t = tagged.begin();
535  typename std::vector< T >::iterator d = data.begin();
536  for( Range::const_iterator i = entities.begin(); i != entities.end() && t != tagged.end(); ++i, d += vals_per_tag )
537  {
538  if( *i == *t )
539  {
540  ++t;
541  rval = mbImpl->tag_get_data( tag, &( *i ), 1, &( *d ) );
542  if( MB_SUCCESS != rval ) return rval;
543  }
544  }
545 
546  // Write the tag values, one entity per line.
547  write_data( stream, data, vals_per_tag );
548 
549  return MB_SUCCESS;
550 }
551 
552 ErrorCode WriteVtk::write_bit_tag( std::ostream& stream, Tag tag, const Range& entities, const Range& tagged )
553 {
554  ErrorCode rval;
555  const unsigned long n = entities.size();
556 
557  // Get tag properties
558 
559  std::string name;
560  int vals_per_tag;
561  if( MB_SUCCESS != mbImpl->tag_get_name( tag, name ) || MB_SUCCESS != mbImpl->tag_get_length( tag, vals_per_tag ) )
562  return MB_FAILURE;
563 
564  if( vals_per_tag > 8 )
565  {
566  MB_SET_ERR( MB_FAILURE, "Invalid tag size for bit tag \"" << name << "\"" );
567  }
568 
569  // Get a tag value for each entity.
570  // Get bits for each entity and unpack into
571  // one integer in the 'data' array for each bit.
572  // Initialize 'data' to zero because we will skip
573  // those entities for which the tag is not set.
574  std::vector< unsigned short > data;
575  data.resize( n * vals_per_tag, 0 );
576  Range::const_iterator t = tagged.begin();
577  std::vector< unsigned short >::iterator d = data.begin();
578  for( Range::const_iterator i = entities.begin(); i != entities.end() && t != tagged.end(); ++i )
579  {
580  if( *i == *t )
581  {
582  ++t;
583  unsigned char value;
584  rval = mbImpl->tag_get_data( tag, &( *i ), 1, &value );
585  for( int j = 0; j < vals_per_tag; ++j, ++d )
586  *d = (unsigned short)( value & ( 1 << j ) ? 1 : 0 );
587  if( MB_SUCCESS != rval ) return rval;
588  }
589  else
590  {
591  // If tag is not set for entity, skip values in array
592  d += vals_per_tag;
593  }
594  }
595 
596  // Write the tag values, one entity per line.
597  write_data( stream, data, vals_per_tag );
598 
599  return MB_SUCCESS;
600 }
601 
602 ErrorCode WriteVtk::write_tag( std::ostream& s, Tag tag, const Range& entities, const Range& tagged )
603 {
604  // Get tag properties
605  std::string name;
606  DataType type;
607  int vals_per_tag;
608  if( MB_SUCCESS != mbImpl->tag_get_name( tag, name ) || MB_SUCCESS != mbImpl->tag_get_length( tag, vals_per_tag ) ||
609  MB_SUCCESS != mbImpl->tag_get_data_type( tag, type ) )
610  return MB_FAILURE;
611 
612  // Skip tags of type ENTITY_HANDLE
613  if( MB_TYPE_HANDLE == type ) return MB_FAILURE;
614 
615  // Now that we're past the point where the name would be used in
616  // an error message, remove any spaces to conform to VTK file.
617  for( std::string::iterator i = name.begin(); i != name.end(); ++i )
618  {
619  if( isspace( *i ) || iscntrl( *i ) ) *i = '_';
620  }
621 
622  // Write the tag description
623  if( 3 == vals_per_tag && MB_TYPE_DOUBLE == type )
624  s << "VECTORS " << name << ' ' << VtkUtil::vtkTypeNames[type] << std::endl;
625  else if( 9 == vals_per_tag )
626  s << "TENSORS " << name << ' ' << VtkUtil::vtkTypeNames[type] << std::endl;
627  else
628  s << "SCALARS " << name << ' ' << VtkUtil::vtkTypeNames[type] << ' ' << vals_per_tag << std::endl
629  << "LOOKUP_TABLE default" << std::endl;
630 
631  // Write the tag data
632  switch( type )
633  {
634  case MB_TYPE_OPAQUE:
635  return write_tag< unsigned char >( s, tag, entities, tagged, 0 );
636  case MB_TYPE_INTEGER:
637  return write_tag< int >( s, tag, entities, tagged, 0 );
638  case MB_TYPE_DOUBLE:
639  return write_tag< double >( s, tag, entities, tagged, 0 );
640  case MB_TYPE_BIT:
641  return write_bit_tag( s, tag, entities, tagged );
642  default:
643  return MB_FAILURE;
644  }
645 }
646 
647 } // namespace moab