moab
moab::WriteSTL Class Reference

ASCII and Binary Stereo Lithography File writers. More...

#include <WriteSTL.hpp>

Inheritance diagram for moab::WriteSTL:
moab::WriterIface

List of all members.

Public Member Functions

 WriteSTL (Interface *impl)
 Constructor.
virtual ~WriteSTL ()
 Destructor.
ErrorCode write_file (const char *file_name, const bool overwrite, const FileOptions &opts, const EntityHandle *output_list, const int num_sets, const std::vector< std::string > &qa_list, const Tag *tag_list=NULL, int num_tags=0, int export_dimension=3)
 writes out a file

Static Public Member Functions

static WriterIfacefactory (Interface *)
 factory method forSTL writer

Protected Types

enum  ByteOrder { STL_BIG_ENDIAN, STL_LITTLE_ENDIAN, STL_UNKNOWN_BYTE_ORDER }

Protected Member Functions

ErrorCode ascii_write_triangles (FILE *file, const char header[82], const Range &triangles, int precision)
 Write list of triangles to an STL file.
ErrorCode binary_write_triangles (FILE *file, const char header[82], ByteOrder byte_order, const Range &triangles)
 Write list of triangles to an STL file.
ErrorCode get_triangle_data (const double vtx_coords[9], float v1[3], float v2[3], float v3[3], float n[3])

Protected Attributes

InterfacembImpl
 interface instance
WriteUtilIfacemWriteIface

Private Member Functions

ErrorCode make_header (char header[82], const std::vector< std::string > &qa_list)
ErrorCode get_triangles (const EntityHandle *set_array, int set_array_length, Range &triangles)
FILE * open_file (const char *name, bool overwrite, bool binary)

Detailed Description

ASCII and Binary Stereo Lithography File writers.

Author:
Jason Kraftcheck

This writer will write only the MBTRI elements in the mesh. It will not decompose other 2-D elements into triangles, nor will it skin the mesh or do any other high-level operation to generate triangles from 3-D elements.

Binary files will be written with a little-endian byte order by default. The byte order can be controlled with writer options.

Definition at line 41 of file WriteSTL.hpp.


Member Enumeration Documentation

enum moab::WriteSTL::ByteOrder [protected]
Enumerator:
STL_BIG_ENDIAN 
STL_LITTLE_ENDIAN 
STL_UNKNOWN_BYTE_ORDER 

Definition at line 68 of file WriteSTL.hpp.


Constructor & Destructor Documentation

Constructor.

Definition at line 57 of file WriteSTL.cpp.

    : mbImpl(impl)
{
  impl->query_interface( mWriteIface );
}
WriteSTL::~WriteSTL ( ) [virtual]

Destructor.

Definition at line 63 of file WriteSTL.cpp.


Member Function Documentation

ErrorCode WriteSTL::ascii_write_triangles ( FILE *  file,
const char  header[82],
const Range triangles,
int  precision 
) [protected]

Write list of triangles to an STL file.

Definition at line 258 of file WriteSTL.cpp.

{
  const char solid_name[] = "MOAB";
  
  char myheader[81] = "solid ";
  strcat( myheader, solid_name );
  strncat( myheader, header, 80 );
  
  if (EOF == fputs( myheader, file ) || EOF == fputs( "\n", file ))
    return MB_FILE_WRITE_ERROR;
  
  ErrorCode rval;
  double coords[9];
  float v1[3], v2[3], v3[3];
  float n[3];
  for (Range::const_iterator iter = triangles.begin();
       iter != triangles.end(); ++iter)
  {
    const EntityHandle* conn;
    int num_vtx;
    
    rval = mbImpl->get_connectivity( *iter, conn, num_vtx );
    if (MB_SUCCESS != rval)
      return rval;
    if (num_vtx != 3)
      return MB_FAILURE;
    
    rval = mbImpl->get_coords( conn, 3, coords );
    if (MB_SUCCESS != rval)
      return rval;
    
    rval = get_triangle_data( coords, v1, v2, v3, n );
    if (MB_SUCCESS != rval)
      return rval;
   
    fprintf( file,"facet normal %e %e %e\n", n[0], n[1], n[2] );
    fprintf( file,"outer loop\n" );
    fprintf( file,"vertex %.*e %.*e %.*e\n", prec, v1[0], prec, v1[1], prec, v1[2] );
    fprintf( file,"vertex %.*e %.*e %.*e\n", prec, v2[0], prec, v2[1], prec, v2[2] );
    fprintf( file,"vertex %.*e %.*e %.*e\n", prec, v3[0], prec, v3[1], prec, v3[2] );
    fprintf( file,"endloop\n" );
    fprintf( file,"endfacet\n" );
  }
  
  fprintf( file,"endsolid %s\n", solid_name );
  return MB_SUCCESS;
}
ErrorCode WriteSTL::binary_write_triangles ( FILE *  file,
const char  header[82],
ByteOrder  byte_order,
const Range triangles 
) [protected]

Write list of triangles to an STL file.

Definition at line 318 of file WriteSTL.cpp.

{
  ErrorCode rval;
  if (fwrite( header, 80, 1, file ) != 1)
    return MB_FILE_WRITE_ERROR;
  
    // default to little endian if byte_order == UNKNOWN_BYTE_ORDER
  const bool want_big_endian = (byte_order == STL_BIG_ENDIAN);
  const bool am_big_endian = !SysUtil::little_endian();
  const bool swap_bytes = (want_big_endian == am_big_endian);
    
  if (triangles.size() > INT_MAX) // can't write that many triangles
    return MB_FAILURE;  
  
  uint32_t count = (uint32_t)triangles.size();
  if (swap_bytes)
    SysUtil::byteswap(&count, 1);
  if (fwrite( &count, 4, 1, file ) != 1)
    return MB_FILE_WRITE_ERROR;

  double coords[9];
  BinTri tri;
  tri.pad[0] = tri.pad[1] = '\0';
  for (Range::const_iterator iter = triangles.begin();
       iter != triangles.end(); ++iter)
  {
    const EntityHandle* conn;
    int num_vtx;
    
    rval = mbImpl->get_connectivity( *iter, conn, num_vtx );
    if (MB_SUCCESS != rval)
      return rval;
    if (num_vtx != 3)
      return MB_FAILURE;
    
    rval = mbImpl->get_coords( conn, 3, coords );
    if (MB_SUCCESS != rval)
      return rval;
    
    rval = get_triangle_data( coords, tri.vertex1, tri.vertex2, tri.vertex3, tri.normal );
    if (MB_SUCCESS != rval)
      return rval;
    
    if (swap_bytes)
    {
      SysUtil::byteswap( tri.normal, 3 );
      SysUtil::byteswap( tri.vertex1, 3 );
      SysUtil::byteswap( tri.vertex2, 3 );
      SysUtil::byteswap( tri.vertex3, 3 );
    }
   
    if (1 != fwrite( &tri, 50, 1, file ))
      return MB_FILE_WRITE_ERROR;
  }
  
  return MB_SUCCESS;
}
WriterIface * WriteSTL::factory ( Interface iface) [static]

factory method forSTL writer

Definition at line 54 of file WriteSTL.cpp.

  { return new WriteSTL( iface ); }
ErrorCode WriteSTL::get_triangle_data ( const double  vtx_coords[9],
float  v1[3],
float  v2[3],
float  v3[3],
float  n[3] 
) [protected]

Given an array of vertex coordinates for a triangle, pass back individual point coordinates as floats and calculate triangle normal.

Definition at line 225 of file WriteSTL.cpp.

{
  float e1[3], e2[3];
  v1[0] = (float)coords[0];
  v1[1] = (float)coords[1];
  v1[2] = (float)coords[2];
  v2[0] = (float)coords[3];
  v2[1] = (float)coords[4];
  v2[2] = (float)coords[5];
  v3[0] = (float)coords[6];
  v3[1] = (float)coords[7];
  v3[2] = (float)coords[8];
  e1[0] = v2[0] - v1[0];
  e1[1] = v2[1] - v1[1];
  e1[2] = v2[2] - v1[2];
  e2[0] = v3[0] - v1[0];
  e2[1] = v3[1] - v1[1];
  e2[2] = v3[2] - v1[2];
  n[0] = e1[1]*e2[2] - e1[2]*e2[1];
  n[1] = e1[2]*e2[0] - e1[0]*e2[2];
  n[2] = e1[0]*e2[1] - e1[1]*e2[0];
  float inv_len = 1.0f / (float)sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
  n[0] *= inv_len;
  n[1] *= inv_len;
  n[2] *= inv_len;
  return MB_SUCCESS;
}
ErrorCode WriteSTL::get_triangles ( const EntityHandle set_array,
int  set_array_length,
Range triangles 
) [private]

Get triangles to write from input array of entity sets. If no sets, gets all triangles.

Definition at line 202 of file WriteSTL.cpp.

{
  if (!set_array || set_array_length == 0)
  {
    return mbImpl->get_entities_by_type( 0, MBTRI, triangles );
  }
  
  const EntityHandle* iter = set_array;
  const EntityHandle* end = iter + set_array_length;
  for (; iter != end; ++iter)
  {
    Range r;
    ErrorCode rval = mbImpl->get_entities_by_type( *iter, MBTRI, r, true );
    if (MB_SUCCESS != rval)
      return rval;
    triangles.merge( r );
  }
  
  return MB_SUCCESS;
}
ErrorCode WriteSTL::make_header ( char  header[82],
const std::vector< std::string > &  qa_list 
) [private]

Construct 80-byte, null-terminated description string from qa_list. Unused space in header will be null-char padded.

Definition at line 182 of file WriteSTL.cpp.

{
  memset( header, 0, 81 );
  
  std::string result;
  for (std::vector<std::string>::const_iterator i = qa_list.begin(); i != qa_list.end(); ++i)
  {
    result += " ";
    result += *i;
  }
  
  size_t len = result.size();
  if (len > 80)
    len = 80;
  memcpy( header, result.c_str(), len );
  
  return MB_SUCCESS;
}
FILE * WriteSTL::open_file ( const char *  name,
bool  overwrite,
bool  binary 
) [private]

Open a file, respecting passed overwrite value and subclass-specified value for need_binary_io().

Definition at line 141 of file WriteSTL.cpp.

{
    // Open file with write access, and create it if it doesn't exist.
  int flags = O_WRONLY|O_CREAT;
    // Select behavior if the named file already exists.  If 
    // overwrite is true, truncate the file.  If it is false,
    // make the call to open() fail.
  if (overwrite)
    flags |= O_TRUNC;
  else
    flags |= O_EXCL;
    // If platform defines a "binary" bit in the file access
    // flags (i.e. we're building on windows), then set it
    // if we're writing a binary file.
#ifdef O_BINARY
  if (binary)
    flags |= O_BINARY;
#endif

    // Give everyone read and write, but not execute, permision.
    // These are not the final permisions for the file.  Permissions
    // are removed according to the user's umask.  All we want to
    // say here is that the executable bits should not be set because
    // this isn't an executable file.  Everything else is a user
    // preference and should be left up to the umask.
  int creat_mode = _S_IREAD|_S_IWRITE;

    // Open the file.
  int fd = open( name, flags, creat_mode );
  if (fd < 0)
  {
    mWriteIface->report_error( "%s: %s\n", name, strerror(errno) );
    return 0;
  }
  FILE* result = fdopen( fd, binary ? "wb": "w" );
  if (!result)
    close( fd );
  
  return result;
}
ErrorCode WriteSTL::write_file ( const char *  file_name,
const bool  overwrite,
const FileOptions opts,
const EntityHandle output_list,
const int  num_sets,
const std::vector< std::string > &  qa_list,
const Tag tag_list = NULL,
int  num_tags = 0,
int  export_dimension = 3 
) [virtual]

writes out a file

Implements moab::WriterIface.

Definition at line 69 of file WriteSTL.cpp.

{
  char header[81];
  Range triangles;
  ErrorCode rval;
  
  if (tag_list && num_tags) {
    mWriteIface->report_error( "STL file does not support tag data\n" );
    return MB_TYPE_OUT_OF_RANGE;
  }
  
  rval = make_header( header, qa_list );
  if (MB_SUCCESS != rval)
    return rval;
  
  rval = get_triangles( ent_handles, num_sets, triangles );
  if (MB_SUCCESS != rval)
    return rval;

  if (triangles.empty()) {
    mWriteIface->report_error( "No triangles to write." );
    return MB_ENTITY_NOT_FOUND;
  }
 
  bool is_ascii = false, is_binary = false;
  if (MB_SUCCESS == opts.get_null_option( "ASCII" ))
    is_ascii = true;
  if (MB_SUCCESS == opts.get_null_option( "BINARY" ))
    is_binary = true;
  if (is_ascii && is_binary) {
    mWriteIface->report_error( "Conflicting options: BINARY ASCII\n" );
    return MB_FAILURE;
  }
  
  bool big_endian = false, little_endian = false;
  if (MB_SUCCESS == opts.get_null_option( "BIG_ENDIAN" ))
    big_endian = true;
  if (MB_SUCCESS == opts.get_null_option( "LITTLE_ENDIAN" ))
    little_endian = true;
  if (big_endian && little_endian) {
    mWriteIface->report_error( "Conflicting options: BIG_ENDIAN LITTLE_ENDIAN\n" );
    return MB_FAILURE;
  }
  ByteOrder byte_order = big_endian ? STL_BIG_ENDIAN : little_endian ? STL_LITTLE_ENDIAN : STL_UNKNOWN_BYTE_ORDER;
    
  FILE* file = open_file( file_name, overwrite, is_binary );
  if (!file)
    return MB_FILE_DOES_NOT_EXIST; 
  
  if (is_binary)
    rval = binary_write_triangles( file, header, byte_order, triangles );
  else {
      // Get precision for node coordinates
    int precision;
    if (MB_SUCCESS != opts.get_int_option( "PRECISION", precision ))
      precision = DEFAULT_PRECISION;

    rval = ascii_write_triangles( file, header, triangles, precision );
  }
  fclose( file );
  return rval;
}

Member Data Documentation

interface instance

Definition at line 91 of file WriteSTL.hpp.

Definition at line 92 of file WriteSTL.hpp.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines