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