libMesh::PostscriptIO Class Reference

#include <postscript_io.h>

Inheritance diagram for libMesh::PostscriptIO:

List of all members.

Public Member Functions

 PostscriptIO (const MeshBase &mesh)
virtual ~PostscriptIO ()
virtual void write (const std::string &)
void plot_quadratic_elem (const Elem *elem)
void plot_linear_elem (const Elem *elem)
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=NULL)
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
unsigned int & ascii_precision ()

Public Attributes

Real shade_value
Real line_width

Protected Member Functions

const MeshBasemesh () const

Protected Attributes

const bool _is_parallel_format

Private Member Functions

void _compute_edge_bezier_coeffs (const Elem *elem)

Private Attributes

std::vector< Point_bezier_coeffs
Point _offset
Real _scale
Point _current_point
std::ostringstream _cell_string
std::ofstream _out

Static Private Attributes

static const float _bezier_transform [3][3]

Detailed Description

This class implements writing 2D meshes in Postscript. It borrows several ideas from, and is a more simple-minded version of, the DataOutBase::write_eps() function from Deal II. Only output is supported here, and only the Mesh (none of the data) is written. The main use I imagined for this class is creating nice Mesh images for publications, since I didn't find/don't know of a free visualization program which would do this.

Author:
John W. Peterson, 2008

Definition at line 51 of file postscript_io.h.


Constructor & Destructor Documentation

libMesh::PostscriptIO::PostscriptIO ( const MeshBase mesh  )  [explicit]

Constructor.

Definition at line 42 of file postscript_io.C.

References _bezier_coeffs.

00043   : MeshOutput<MeshBase> (mesh_in),
00044     shade_value(0.0),
00045     line_width(0.5),
00046     //_M(3,3),
00047     _offset(0., 0.),
00048     _scale(1.0),
00049     _current_point(0., 0.)
00050 {
00051   // This code is still undergoing some development.
00052   libmesh_experimental();
00053 
00054   // Entries of transformation matrix from physical to Bezier coords.
00055   // _M(0,0) = -1./6.;    _M(0,1) = 1./6.;    _M(0,2) = 1.;
00056   // _M(1,0) = -1./6.;    _M(1,1) = 0.5  ;    _M(1,2) = 1./6.;
00057   // _M(2,0) = 0.    ;    _M(2,1) = 1.   ;    _M(2,2) = 0.;
00058 
00059   // Make sure there is enough room to store Bezier coefficients.
00060   _bezier_coeffs.resize(3);
00061 }

libMesh::PostscriptIO::~PostscriptIO (  )  [virtual]

Destructor.

Definition at line 65 of file postscript_io.C.

00066 {
00067 }


Member Function Documentation

void libMesh::PostscriptIO::_compute_edge_bezier_coeffs ( const Elem elem  )  [private]

Given a quadratic edge Elem which lies in the x-y plane, computes the Bezier coefficients. These may be passed to the Postscript routine "curveto".

Definition at line 271 of file postscript_io.C.

References _bezier_coeffs, _bezier_transform, _offset, _scale, libMeshEnums::EDGE3, libMesh::Elem::point(), and libMesh::Elem::type().

Referenced by plot_quadratic_elem().

00272 {
00273   // I only know how to do this for an Edge3!
00274   libmesh_assert_equal_to (elem->type(), EDGE3);
00275 
00276   // Get x-coordinates into an array, transform them,
00277   // and repeat for y.
00278   float
00279     phys_coords[3] = {0., 0., 0.},
00280     bez_coords[3]  = {0., 0., 0.};
00281 
00282   for (unsigned int i=0; i<2; ++i)
00283     {
00284       // Initialize vectors.  Physical coordinates are initialized
00285       // by their postscript-scaled values.
00286       for (unsigned int j=0; j<3; ++j)
00287         {
00288           phys_coords[j] = (elem->point(j)(i) - _offset(i)) * _scale;
00289           bez_coords[j] = 0.; // zero out result vector
00290         }
00291 
00292       // Multiply matrix times vector
00293       for (unsigned int j=0; j<3; ++j)
00294         for (unsigned int k=0; k<3; ++k)
00295           bez_coords[j] += _bezier_transform[j][k]*phys_coords[k];
00296 
00297       // Store result in _bezier_coeffs
00298       for (unsigned int j=0; j<3; ++j)
00299         _bezier_coeffs[j](i) = phys_coords[j];
00300     }
00301 }

unsigned int& libMesh::MeshOutput< MeshBase >::ascii_precision (  )  [inherited]

Return/set the precision to use when writing ASCII files.

By default we use numeric_limits<Real>::digits10 + 2, which should be enough to write out to ASCII and get the exact same Real back when reading in.

Referenced by libMesh::TecplotIO::write_ascii(), libMesh::GMVIO::write_ascii_new_impl(), and libMesh::GMVIO::write_ascii_old_impl().

void libMesh::PostscriptIO::plot_linear_elem ( const Elem elem  ) 

Draws an element with straight lines

Definition at line 197 of file postscript_io.C.

References _cell_string, _current_point, _offset, _out, _scale, libMesh::Elem::n_vertices(), libMesh::Elem::point(), and shade_value.

Referenced by write().

00198 {
00199   // Clear the string contents.  Yes, this really is how you do that...
00200   _cell_string.str("");
00201 
00202   // The general strategy is:
00203   // 1.) Use m  := {moveto} to go to vertex 0.
00204   // 2.) Use l  := {lineto} commands to draw lines to vertex 1, 2, ... N-1.
00205   // 3a.) Use lx := {lineto closepath stroke} command at  vertex N to draw the last line.
00206   // 3b.)     lf := {lineto closepath fill} command to shade the cell just drawn
00207   // All of our 2D elements' vertices are numbered in counterclockwise order,
00208   // so we can just draw them in the same order.
00209 
00210   // 1.)
00211   _current_point = (elem->point(0) - _offset) * _scale;
00212   _cell_string << _current_point(0) << " " << _current_point(1) << " "; // write x y
00213   _cell_string << "m ";
00214 
00215   // 2.)
00216   const unsigned int nv=elem->n_vertices();
00217   for (unsigned int v=1; v<nv-1; ++v)
00218     {
00219       _current_point = (elem->point(v) - _offset) * _scale;
00220       _cell_string << _current_point(0) << " " << _current_point(1) << " "; // write x y
00221       _cell_string << "l ";
00222     }
00223 
00224   // 3.)
00225   _current_point = (elem->point(nv-1) - _offset) * _scale;
00226   _cell_string << _current_point(0) << " " << _current_point(1) << " "; // write x y
00227 
00228   // We draw the shaded (interior) parts first, if applicable.
00229   if (shade_value > 0.0)
00230     _out << shade_value << " sg " << _cell_string.str() << "lf\n";
00231 
00232   // Draw the black lines (I guess we will always do this)
00233   _out << "0 sg " << _cell_string.str() << "lx\n";
00234 }

void libMesh::PostscriptIO::plot_quadratic_elem ( const Elem elem  ) 

Draws an element with Bezier curves

Definition at line 239 of file postscript_io.C.

References _bezier_coeffs, _compute_edge_bezier_coeffs(), _current_point, _offset, _out, _scale, libMesh::Elem::build_side(), libMeshEnums::EDGE3, libMesh::AutoPtr< Tp >::get(), libMesh::Elem::n_sides(), and side.

00240 {
00241   for (unsigned int ns=0; ns<elem->n_sides(); ++ns)
00242     {
00243       // Build the quadratic side
00244       AutoPtr<Elem> side = elem->build_side(ns);
00245 
00246       // Be sure it's quadratic (Edge2).  Eventually we could
00247       // handle cubic elements as well...
00248       libmesh_assert_equal_to ( side->type(), EDGE3 );
00249 
00250       _out << "0 sg ";
00251 
00252       // Move to the first point on this side.
00253       _current_point = (side->point(0) - _offset) * _scale;
00254       _out << _current_point(0) << " " << _current_point(1) << " "; // write x y
00255       _out << "m ";
00256 
00257       // Compute _bezier_coeffs for this edge.  This fills up
00258       // the _bezier_coeffs vector.
00259       this->_compute_edge_bezier_coeffs(side.get());
00260 
00261       // Print curveto path to file
00262       for (unsigned int i=0; i<_bezier_coeffs.size(); ++i)
00263         _out << _bezier_coeffs[i](0) << " " << _bezier_coeffs[i](1) << " ";
00264       _out << " cs\n";
00265     }
00266 }

void libMesh::PostscriptIO::write ( const std::string &  fname  )  [virtual]

This method implements writing a mesh to a specified file.

Implements libMesh::MeshOutput< MeshBase >.

Definition at line 71 of file postscript_io.C.

References libMesh::MeshOutput< MeshBase >::_is_parallel_format, _offset, _out, _scale, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshTools::bounding_box(), line_width, libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), plot_linear_elem(), libMesh::processor_id(), and libMesh::Real.

00072 {
00073   // We may need to gather a ParallelMesh to output it, making that
00074   // const qualifier in our constructor a dirty lie
00075   MeshSerializer serialize(const_cast<MeshBase&>(this->mesh()), !_is_parallel_format);
00076 
00077   if (libMesh::processor_id() == 0)
00078     {
00079       // Get a constant reference to the mesh.
00080       const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00081 
00082       // Only works in 2D
00083       libmesh_assert_equal_to (the_mesh.mesh_dimension(), 2);
00084 
00085       // Create output file stream.
00086       // _out is now a private member of the class.
00087       _out.open(fname.c_str());
00088 
00089       // Make sure it opened correctly
00090       if (!_out.good())
00091         libmesh_file_error(fname.c_str());
00092 
00093       // The mesh bounding box gives us info about what the
00094       // Postscript bounding box should be.
00095       MeshTools::BoundingBox bbox = MeshTools::bounding_box(the_mesh);
00096 
00097       // Add a little extra padding to the "true" bounding box so
00098       // that we can still see the boundary
00099       const Real percent_padding = 0.01;
00100       const Real dx=bbox.second(0)-bbox.first(0); libmesh_assert_greater (dx, 0.0);
00101       const Real dy=bbox.second(1)-bbox.first(1); libmesh_assert_greater (dy, 0.0);
00102 
00103       const Real x_min = bbox.first(0)  - percent_padding*dx;
00104       const Real y_min = bbox.first(1)  - percent_padding*dy;
00105       const Real x_max = bbox.second(0) + percent_padding*dx;
00106       const Real y_max = bbox.second(1) + percent_padding*dy;
00107 
00108       // Width of the output as given in postscript units.
00109       // This usually is given by the strange unit 1/72 inch.
00110       // A width of 300 represents a size of roughly 10 cm.
00111       const Real width = 300;
00112       _scale = width / (x_max-x_min);
00113       _offset(0) = x_min;
00114       _offset(1) = y_min;
00115 
00116       // Header writing stuff stolen from Deal.II
00117       std::time_t  time1= std::time (0);
00118       std::tm     *time = std::localtime(&time1);
00119       _out << "%!PS-Adobe-2.0 EPSF-1.2" << '\n'
00120         //<< "%!PS-Adobe-1.0" << '\n' // Lars' PS version
00121           << "%%Filename: " << fname << '\n'
00122           << "%%Title: LibMesh Output" << '\n'
00123           << "%%Creator: LibMesh: A C++ finite element library" << '\n'
00124           << "%%Creation Date: "
00125           << time->tm_year+1900 << "/"
00126           << time->tm_mon+1 << "/"
00127           << time->tm_mday << " - "
00128           << time->tm_hour << ":"
00129           << std::setw(2) << time->tm_min << ":"
00130           << std::setw(2) << time->tm_sec << '\n'
00131           << "%%BoundingBox: "
00132         // lower left corner
00133           << "0 0 "
00134         // upper right corner
00135           << static_cast<unsigned int>( rint((x_max-x_min) * _scale ))
00136           << ' '
00137           << static_cast<unsigned int>( rint((y_max-y_min) * _scale ))
00138           << '\n';
00139 
00140       // define some abbreviations to keep
00141       // the output small:
00142       // m=move turtle to
00143       // l=define a line
00144       // s=set rgb color
00145       // sg=set gray value
00146       // lx=close the line and plot the line
00147       // lf=close the line and fill the interior
00148       _out << "/m {moveto} bind def"      << '\n'
00149           << "/l {lineto} bind def"      << '\n'
00150           << "/s {setrgbcolor} bind def" << '\n'
00151           << "/sg {setgray} bind def"    << '\n'
00152           << "/cs {curveto stroke} bind def" << '\n'
00153           << "/lx {lineto closepath stroke} bind def" << '\n'
00154           << "/lf {lineto closepath fill} bind def"   << '\n';
00155 
00156       _out << "%%EndProlog" << '\n';
00157       //          << '\n';
00158 
00159       // Set line width in the postscript file.
00160       _out << line_width << " setlinewidth" << '\n';
00161 
00162       // Set line cap and join options
00163       _out << "1 setlinecap" << '\n';
00164       _out << "1 setlinejoin" << '\n';
00165 
00166       // allow only five digits for output (instead of the default
00167       // six); this should suffice even for fine grids, but reduces
00168       // the file size significantly
00169       _out << std::setprecision (5);
00170 
00171       // Loop over the active elements, draw lines for the edges.  We
00172       // draw even quadratic elements with straight sides, i.e. a straight
00173       // line sits between each pair of vertices.  Also we draw every edge
00174       // for an element regardless of the fact that it may overlap with
00175       // another.  This would probably be a useful optimization...
00176       MeshBase::const_element_iterator       el     = the_mesh.active_elements_begin();
00177       const MeshBase::const_element_iterator end_el = the_mesh.active_elements_end();
00178       for ( ; el != end_el; ++el)
00179         {
00180           //const Elem* elem = *el;
00181 
00182           this->plot_linear_elem(*el);
00183           //this->plot_quadratic_elem(*el); // Experimental
00184         }
00185 
00186       // Issue the showpage command, and we're done.
00187       _out << "showpage" << std::endl;
00188 
00189     } // end if (libMesh::processor_id() == 0)
00190 }

virtual void libMesh::MeshOutput< MeshBase >::write_equation_systems ( const std::string &  ,
const EquationSystems ,
const std::set< std::string > *  system_names = NULL 
) [virtual, inherited]

This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object.

Referenced by libMesh::Nemesis_IO::write_timestep().

virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  ,
const std::vector< Number > &  ,
const std::vector< std::string > &   
) [inline, virtual, inherited]

This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided.

Reimplemented in libMesh::ExodusII_IO, libMesh::GmshIO, libMesh::GMVIO, libMesh::GnuPlotIO, libMesh::MEDITIO, libMesh::Nemesis_IO, libMesh::TecplotIO, libMesh::UCDIO, and libMesh::VTKIO.

Definition at line 98 of file mesh_output.h.

00101   { libmesh_error(); }


Member Data Documentation

Vector containing 3 points corresponding to Bezier coefficients, as computed by _compute_edge_bezier_coeffs.

Definition at line 118 of file postscript_io.h.

Referenced by _compute_edge_bezier_coeffs(), plot_quadratic_elem(), and PostscriptIO().

const float libMesh::PostscriptIO::_bezier_transform [static, private]
Initial value:
{
  {-1./6., 1./6.,    1.},
  {-1./6.,   0.5, 1./6.},
  {    0.,    1.,    0.}
}

Coefficients of the transformation from physical-space edge coordinates to Bezier basis coefficients. Transforms x and y separately.

Definition at line 112 of file postscript_io.h.

Referenced by _compute_edge_bezier_coeffs().

std::ostringstream libMesh::PostscriptIO::_cell_string [private]

Drawing style-independent data for a single cell. This can be used as a temporary buffer for storing data which may be sent to the output stream multiple times.

Definition at line 140 of file postscript_io.h.

Referenced by plot_linear_elem().

A point object used for temporary calculations

Definition at line 133 of file postscript_io.h.

Referenced by plot_linear_elem(), and plot_quadratic_elem().

const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format [protected, inherited]

Flag specifying whether this format is parallel-capable. If this is false (default) I/O is only permitted when the mesh has been serialized.

Definition at line 126 of file mesh_output.h.

Referenced by write(), libMesh::FroIO::write(), libMesh::EnsightIO::write(), and libMesh::DivaIO::write().

Amount to add to every (x,y) point to place it in Postscript coordinates.

Definition at line 123 of file postscript_io.h.

Referenced by _compute_edge_bezier_coeffs(), plot_linear_elem(), plot_quadratic_elem(), and write().

std::ofstream libMesh::PostscriptIO::_out [private]

Output file stream which will be opened when the file name is known

Definition at line 145 of file postscript_io.h.

Referenced by plot_linear_elem(), plot_quadratic_elem(), and write().

Amount by which to stretch each point to place it in Postscript coordinates.

Definition at line 128 of file postscript_io.h.

Referenced by _compute_edge_bezier_coeffs(), plot_linear_elem(), plot_quadratic_elem(), and write().

Control the thickness of the lines used. 0.5 is a reasonable default for printed images, but you may need to decrease this value (or choose it adaptively) when there are very slim cells present in the mesh.

Definition at line 86 of file postscript_io.h.

Referenced by write().

Controls greyscale shading of cells. By default this value is 0.0 (which actually corresponds to black) and this indicates "no shading" i.e. only mesh lines will be drawn. Any other value in (0,1] will cause the cells to be grey-shaded to some degree, with higher values being lighter. A value of 0.75 gives decent results.

Definition at line 78 of file postscript_io.h.

Referenced by plot_linear_elem().


The documentation for this class was generated from the following files:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:34 UTC

Hosted By:
SourceForge.net Logo