libMesh::PostscriptIO Class Reference

#include <postscript_io.h>

Inheritance diagram for libMesh::PostscriptIO:

Public Member Functions

 PostscriptIO (const MeshBase &mesh)
 
virtual ~PostscriptIO ()
 
virtual void write (const std::string &) override
 
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=nullptr)
 
virtual void write_discontinuous_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
 
virtual void write_nodal_data (const std::string &, const NumericVector< Number > &, const std::vector< std::string > &)
 
virtual void write_nodal_data_discontinuous (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
 
const bool _serial_only_needed_on_proc_0
 

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
Date
2008

Definition at line 53 of file postscript_io.h.

Constructor & Destructor Documentation

◆ PostscriptIO()

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

Constructor.

Definition at line 42 of file postscript_io.C.

References _bezier_coeffs.

42  :
43  MeshOutput<MeshBase> (mesh_in),
44  shade_value(0.0),
45  line_width(0.5),
46  //_M(3,3),
47  _offset(0., 0.),
48  _scale(1.0),
49  _current_point(0., 0.)
50 {
51  // This code is still undergoing some development.
52  libmesh_experimental();
53 
54  // Entries of transformation matrix from physical to Bezier coords.
55  // _M(0,0) = -1./6.; _M(0,1) = 1./6.; _M(0,2) = 1.;
56  // _M(1,0) = -1./6.; _M(1,1) = 0.5 ; _M(1,2) = 1./6.;
57  // _M(2,0) = 0. ; _M(2,1) = 1. ; _M(2,2) = 0.;
58 
59  // Make sure there is enough room to store Bezier coefficients.
60  _bezier_coeffs.resize(3);
61 }
std::vector< Point > _bezier_coeffs

◆ ~PostscriptIO()

libMesh::PostscriptIO::~PostscriptIO ( )
virtual

Destructor.

Definition at line 65 of file postscript_io.C.

66 {
67 }

Member Function Documentation

◆ _compute_edge_bezier_coeffs()

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 264 of file postscript_io.C.

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

Referenced by plot_quadratic_elem().

265 {
266  // I only know how to do this for an Edge3!
267  libmesh_assert_equal_to (elem->type(), EDGE3);
268 
269  // Get x-coordinates into an array, transform them,
270  // and repeat for y.
271  float phys_coords[3] = {0., 0., 0.};
272  float bez_coords[3] = {0., 0., 0.};
273 
274  for (unsigned int i=0; i<2; ++i)
275  {
276  // Initialize vectors. Physical coordinates are initialized
277  // by their postscript-scaled values.
278  for (unsigned int j=0; j<3; ++j)
279  {
280  phys_coords[j] = static_cast<float>
281  ((elem->point(j)(i) - _offset(i)) * _scale);
282  bez_coords[j] = 0.; // zero out result vector
283  }
284 
285  // Multiply matrix times vector
286  for (unsigned int j=0; j<3; ++j)
287  for (unsigned int k=0; k<3; ++k)
288  bez_coords[j] += _bezier_transform[j][k]*phys_coords[k];
289 
290  // Store result in _bezier_coeffs
291  for (unsigned int j=0; j<3; ++j)
292  _bezier_coeffs[j](i) = phys_coords[j];
293  }
294 }
static const float _bezier_transform[3][3]
std::vector< Point > _bezier_coeffs

◆ ascii_precision()

unsigned int & libMesh::MeshOutput< MeshBase >::ascii_precision ( )
inlineinherited

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.

Definition at line 244 of file mesh_output.h.

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

245 {
246  return _ascii_precision;
247 }

◆ mesh()

◆ plot_linear_elem()

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

Draws an element with straight lines

Definition at line 190 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().

191 {
192  // Clear the string contents. Yes, this really is how you do that...
193  _cell_string.str("");
194 
195  // The general strategy is:
196  // 1.) Use m := {moveto} to go to vertex 0.
197  // 2.) Use l := {lineto} commands to draw lines to vertex 1, 2, ... N-1.
198  // 3a.) Use lx := {lineto closepath stroke} command at vertex N to draw the last line.
199  // 3b.) lf := {lineto closepath fill} command to shade the cell just drawn
200  // All of our 2D elements' vertices are numbered in counterclockwise order,
201  // so we can just draw them in the same order.
202 
203  // 1.)
204  _current_point = (elem->point(0) - _offset) * _scale;
205  _cell_string << _current_point(0) << " " << _current_point(1) << " "; // write x y
206  _cell_string << "m ";
207 
208  // 2.)
209  const unsigned int nv=elem->n_vertices();
210  for (unsigned int v=1; v<nv-1; ++v)
211  {
212  _current_point = (elem->point(v) - _offset) * _scale;
213  _cell_string << _current_point(0) << " " << _current_point(1) << " "; // write x y
214  _cell_string << "l ";
215  }
216 
217  // 3.)
218  _current_point = (elem->point(nv-1) - _offset) * _scale;
219  _cell_string << _current_point(0) << " " << _current_point(1) << " "; // write x y
220 
221  // We draw the shaded (interior) parts first, if applicable.
222  if (shade_value > 0.0)
223  _out << shade_value << " sg " << _cell_string.str() << "lf\n";
224 
225  // Draw the black lines (I guess we will always do this)
226  _out << "0 sg " << _cell_string.str() << "lx\n";
227 }
std::ostringstream _cell_string

◆ plot_quadratic_elem()

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

Draws an element with Bezier curves

Definition at line 232 of file postscript_io.C.

References _bezier_coeffs, _compute_edge_bezier_coeffs(), _current_point, _offset, _out, _scale, libMesh::Elem::build_side_ptr(), libMesh::EDGE3, side, and libMesh::Elem::side_index_range().

233 {
234  for (auto ns : elem->side_index_range())
235  {
236  // Build the quadratic side
237  std::unique_ptr<const Elem> side = elem->build_side_ptr(ns);
238 
239  // Be sure it's quadratic (Edge2). Eventually we could
240  // handle cubic elements as well...
241  libmesh_assert_equal_to ( side->type(), EDGE3 );
242 
243  _out << "0 sg ";
244 
245  // Move to the first point on this side.
246  _current_point = (side->point(0) - _offset) * _scale;
247  _out << _current_point(0) << " " << _current_point(1) << " "; // write x y
248  _out << "m ";
249 
250  // Compute _bezier_coeffs for this edge. This fills up
251  // the _bezier_coeffs vector.
252  this->_compute_edge_bezier_coeffs(side.get());
253 
254  // Print curveto path to file
255  for (std::size_t i=0; i<_bezier_coeffs.size(); ++i)
256  _out << _bezier_coeffs[i](0) << " " << _bezier_coeffs[i](1) << " ";
257  _out << " cs\n";
258  }
259 }
void _compute_edge_bezier_coeffs(const Elem *elem)
unsigned short int side
Definition: xdr_io.C:50
std::vector< Point > _bezier_coeffs

◆ write()

void libMesh::PostscriptIO::write ( const std::string &  fname)
overridevirtual

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_element_ptr_range(), libMesh::MeshTools::create_bounding_box(), line_width, libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::mesh_dimension(), plot_linear_elem(), and libMesh::Real.

72 {
73  // We may need to gather a DistributedMesh to output it, making that
74  // const qualifier in our constructor a dirty lie
75  MeshSerializer serialize(const_cast<MeshBase &>(this->mesh()), !_is_parallel_format);
76 
77  if (this->mesh().processor_id() == 0)
78  {
79  // Get a constant reference to the mesh.
80  const MeshBase & the_mesh = MeshOutput<MeshBase>::mesh();
81 
82  // Only works in 2D
83  libmesh_assert_equal_to (the_mesh.mesh_dimension(), 2);
84 
85  // Create output file stream.
86  // _out is now a private member of the class.
87  _out.open(fname.c_str());
88 
89  // Make sure it opened correctly
90  if (!_out.good())
91  libmesh_file_error(fname.c_str());
92 
93  // The mesh bounding box gives us info about what the
94  // Postscript bounding box should be.
95  BoundingBox bbox = MeshTools::create_bounding_box(the_mesh);
96 
97  // Add a little extra padding to the "true" bounding box so
98  // that we can still see the boundary
99  const Real percent_padding = 0.01;
100  const Real dx=bbox.second(0)-bbox.first(0); libmesh_assert_greater (dx, 0.0);
101  const Real dy=bbox.second(1)-bbox.first(1); libmesh_assert_greater (dy, 0.0);
102 
103  const Real x_min = bbox.first(0) - percent_padding*dx;
104  const Real y_min = bbox.first(1) - percent_padding*dy;
105  const Real x_max = bbox.second(0) + percent_padding*dx;
106  const Real y_max = bbox.second(1) + percent_padding*dy;
107 
108  // Width of the output as given in postscript units.
109  // This usually is given by the strange unit 1/72 inch.
110  // A width of 300 represents a size of roughly 10 cm.
111  const Real width = 300;
112  _scale = width / (x_max-x_min);
113  _offset(0) = x_min;
114  _offset(1) = y_min;
115 
116  // Header writing stuff stolen from Deal.II
117  std::time_t time1= std::time (0);
118  std::tm * time = std::localtime(&time1);
119  _out << "%!PS-Adobe-2.0 EPSF-1.2" << '\n'
120  //<< "%!PS-Adobe-1.0" << '\n' // Lars' PS version
121  << "%%Filename: " << fname << '\n'
122  << "%%Title: LibMesh Output" << '\n'
123  << "%%Creator: LibMesh: A C++ finite element library" << '\n'
124  << "%%Creation Date: "
125  << time->tm_year+1900 << "/"
126  << time->tm_mon+1 << "/"
127  << time->tm_mday << " - "
128  << time->tm_hour << ":"
129  << std::setw(2) << time->tm_min << ":"
130  << std::setw(2) << time->tm_sec << '\n'
131  << "%%BoundingBox: "
132  // lower left corner
133  << "0 0 "
134  // upper right corner
135  << static_cast<unsigned int>( rint((x_max-x_min) * _scale ))
136  << ' '
137  << static_cast<unsigned int>( rint((y_max-y_min) * _scale ))
138  << '\n';
139 
140  // define some abbreviations to keep
141  // the output small:
142  // m=move turtle to
143  // l=define a line
144  // s=set rgb color
145  // sg=set gray value
146  // lx=close the line and plot the line
147  // lf=close the line and fill the interior
148  _out << "/m {moveto} bind def" << '\n'
149  << "/l {lineto} bind def" << '\n'
150  << "/s {setrgbcolor} bind def" << '\n'
151  << "/sg {setgray} bind def" << '\n'
152  << "/cs {curveto stroke} bind def" << '\n'
153  << "/lx {lineto closepath stroke} bind def" << '\n'
154  << "/lf {lineto closepath fill} bind def" << '\n';
155 
156  _out << "%%EndProlog" << '\n';
157  // << '\n';
158 
159  // Set line width in the postscript file.
160  _out << line_width << " setlinewidth" << '\n';
161 
162  // Set line cap and join options
163  _out << "1 setlinecap" << '\n';
164  _out << "1 setlinejoin" << '\n';
165 
166  // allow only five digits for output (instead of the default
167  // six); this should suffice even for fine grids, but reduces
168  // the file size significantly
169  _out << std::setprecision (5);
170 
171  // Loop over the active elements, draw lines for the edges. We
172  // draw even quadratic elements with straight sides, i.e. a straight
173  // line sits between each pair of vertices. Also we draw every edge
174  // for an element regardless of the fact that it may overlap with
175  // another. This would probably be a useful optimization...
176  for (const auto & elem : the_mesh.active_element_ptr_range())
177  this->plot_linear_elem(elem);
178 
179  // Issue the showpage command, and we're done.
180  _out << "showpage" << std::endl;
181 
182  } // end if (this->mesh().processor_id() == 0)
183 }
const MeshBase & mesh() const
Definition: mesh_output.h:234
void plot_linear_elem(const Elem *elem)
libMesh::BoundingBox create_bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:386
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ write_discontinuous_equation_systems()

void libMesh::MeshOutput< MeshBase >::write_discontinuous_equation_systems ( const std::string &  fname,
const EquationSystems es,
const std::set< std::string > *  system_names = nullptr 
)
virtualinherited

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

Definition at line 92 of file mesh_output.C.

References libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), and libMesh::out.

Referenced by libMesh::ExodusII_IO::write_timestep_discontinuous().

95 {
96  LOG_SCOPE("write_discontinuous_equation_systems()", "MeshOutput");
97 
98  // We may need to gather and/or renumber a DistributedMesh to output
99  // it, making that const qualifier in our constructor a dirty lie
100  MT & my_mesh = const_cast<MT &>(*_obj);
101 
102  // If we're asked to write data that's associated with a different
103  // mesh, output files full of garbage are the result.
104  libmesh_assert_equal_to(&es.get_mesh(), _obj);
105 
106  // A non-renumbered mesh may not have a contiguous numbering, and
107  // that needs to be fixed before we can build a solution vector.
108  if (my_mesh.max_elem_id() != my_mesh.n_elem() ||
109  my_mesh.max_node_id() != my_mesh.n_nodes())
110  {
111  // If we were allowed to renumber then we should have already
112  // been properly renumbered...
113  libmesh_assert(!my_mesh.allow_renumbering());
114 
115  libmesh_do_once(libMesh::out <<
116  "Warning: This MeshOutput subclass only supports meshes which are contiguously renumbered!"
117  << std::endl;);
118 
119  my_mesh.allow_renumbering(true);
120 
121  my_mesh.renumber_nodes_and_elements();
122 
123  // Not sure what good going back to false will do here, the
124  // renumbering horses have already left the barn...
125  my_mesh.allow_renumbering(false);
126  }
127 
128  MeshSerializer serialize(const_cast<MT &>(*_obj), !_is_parallel_format, _serial_only_needed_on_proc_0);
129 
130  // Build the list of variable names that will be written.
131  std::vector<std::string> names;
132  es.build_variable_names (names, nullptr, system_names);
133 
134  if (!_is_parallel_format)
135  {
136  // Build the nodal solution values & get the variable
137  // names from the EquationSystems object
138  std::vector<Number> soln;
139  es.build_discontinuous_solution_vector (soln, system_names);
140 
141  this->write_nodal_data_discontinuous (fname, soln, names);
142  }
143  else // _is_parallel_format
144  {
145  libmesh_not_implemented();
146  }
147 }
virtual void write_nodal_data_discontinuous(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
Definition: mesh_output.h:114
const MeshBase *const _obj
Definition: mesh_output.h:177
OStreamProxy out(std::cout)

◆ write_equation_systems()

void libMesh::MeshOutput< MeshBase >::write_equation_systems ( const std::string &  fname,
const EquationSystems es,
const std::set< std::string > *  system_names = nullptr 
)
virtualinherited

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

Reimplemented in libMesh::NameBasedIO.

Definition at line 31 of file mesh_output.C.

References libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::EquationSystems::build_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), and libMesh::out.

Referenced by libMesh::Nemesis_IO::write_timestep(), and libMesh::ExodusII_IO::write_timestep().

34 {
35  LOG_SCOPE("write_equation_systems()", "MeshOutput");
36 
37  // We may need to gather and/or renumber a DistributedMesh to output
38  // it, making that const qualifier in our constructor a dirty lie
39  MT & my_mesh = const_cast<MT &>(*_obj);
40 
41  // If we're asked to write data that's associated with a different
42  // mesh, output files full of garbage are the result.
43  libmesh_assert_equal_to(&es.get_mesh(), _obj);
44 
45  // A non-renumbered mesh may not have a contiguous numbering, and
46  // that needs to be fixed before we can build a solution vector.
47  if (my_mesh.max_elem_id() != my_mesh.n_elem() ||
48  my_mesh.max_node_id() != my_mesh.n_nodes())
49  {
50  // If we were allowed to renumber then we should have already
51  // been properly renumbered...
52  libmesh_assert(!my_mesh.allow_renumbering());
53 
54  libmesh_do_once(libMesh::out <<
55  "Warning: This MeshOutput subclass only supports meshes which are contiguously renumbered!"
56  << std::endl;);
57 
58  my_mesh.allow_renumbering(true);
59 
60  my_mesh.renumber_nodes_and_elements();
61 
62  // Not sure what good going back to false will do here, the
63  // renumbering horses have already left the barn...
64  my_mesh.allow_renumbering(false);
65  }
66 
67  MeshSerializer serialize(const_cast<MT &>(*_obj), !_is_parallel_format, _serial_only_needed_on_proc_0);
68 
69  // Build the list of variable names that will be written.
70  std::vector<std::string> names;
71  es.build_variable_names (names, nullptr, system_names);
72 
74  {
75  // Build the nodal solution values & get the variable
76  // names from the EquationSystems object
77  std::vector<Number> soln;
78  es.build_solution_vector (soln, system_names);
79 
80  this->write_nodal_data (fname, soln, names);
81  }
82  else // _is_parallel_format
83  {
84  std::unique_ptr<NumericVector<Number>> parallel_soln =
85  es.build_parallel_solution_vector(system_names);
86 
87  this->write_nodal_data (fname, *parallel_soln, names);
88  }
89 }
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
Definition: mesh_output.h:105
const MeshBase *const _obj
Definition: mesh_output.h:177
OStreamProxy out(std::cout)

◆ write_nodal_data() [1/2]

virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  ,
const std::vector< Number > &  ,
const std::vector< std::string > &   
)
inlinevirtualinherited

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::Nemesis_IO, libMesh::UCDIO, libMesh::NameBasedIO, libMesh::GmshIO, libMesh::GMVIO, libMesh::VTKIO, libMesh::MEDITIO, libMesh::GnuPlotIO, and libMesh::TecplotIO.

Definition at line 105 of file mesh_output.h.

108  { libmesh_not_implemented(); }

◆ write_nodal_data() [2/2]

void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  fname,
const NumericVector< Number > &  parallel_soln,
const std::vector< std::string > &  names 
)
virtualinherited

This method should be overridden by "parallel" output formats for writing nodal data. Instead of getting a localized copy of the nodal solution vector, it is passed a NumericVector of type=PARALLEL which is in node-major order i.e. (u0,v0,w0, u1,v1,w1, u2,v2,w2, u3,v3,w3, ...) and contains n_nodes*n_vars total entries. Then, it is up to the individual I/O class to extract the required solution values from this vector and write them in parallel.

If not implemented, localizes the parallel vector into a std::vector and calls the other version of this function.

Reimplemented in libMesh::Nemesis_IO.

Definition at line 150 of file mesh_output.C.

References libMesh::NumericVector< T >::localize().

153 {
154  // This is the fallback implementation for parallel I/O formats that
155  // do not yet implement proper writing in parallel, and instead rely
156  // on the full solution vector being available on all processors.
157  std::vector<Number> soln;
158  parallel_soln.localize(soln);
159  this->write_nodal_data(fname, soln, names);
160 }
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
Definition: mesh_output.h:105
virtual void localize(std::vector< T > &v_local) const =0

◆ write_nodal_data_discontinuous()

virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data_discontinuous ( const std::string &  ,
const std::vector< Number > &  ,
const std::vector< std::string > &   
)
inlinevirtualinherited

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

Reimplemented in libMesh::ExodusII_IO.

Definition at line 114 of file mesh_output.h.

117  { libmesh_not_implemented(); }

Member Data Documentation

◆ _bezier_coeffs

std::vector<Point> libMesh::PostscriptIO::_bezier_coeffs
private

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

Definition at line 120 of file postscript_io.h.

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

◆ _bezier_transform

const float libMesh::PostscriptIO::_bezier_transform
staticprivate
Initial value:
=
{
{-1.f/6.f, 1.f/6.f, 1.},
{-1.f/6.f, 0.5, 1.f/6.f},
{0., 1., 0.}
}

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

Definition at line 114 of file postscript_io.h.

Referenced by _compute_edge_bezier_coeffs().

◆ _cell_string

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 142 of file postscript_io.h.

Referenced by plot_linear_elem().

◆ _current_point

Point libMesh::PostscriptIO::_current_point
private

A point object used for temporary calculations

Definition at line 135 of file postscript_io.h.

Referenced by plot_linear_elem(), and plot_quadratic_elem().

◆ _is_parallel_format

const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format
protectedinherited

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 159 of file mesh_output.h.

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

◆ _offset

Point libMesh::PostscriptIO::_offset
private

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

Definition at line 125 of file postscript_io.h.

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

◆ _out

std::ofstream libMesh::PostscriptIO::_out
private

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

Definition at line 147 of file postscript_io.h.

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

◆ _scale

Real libMesh::PostscriptIO::_scale
private

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

Definition at line 130 of file postscript_io.h.

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

◆ _serial_only_needed_on_proc_0

const bool libMesh::MeshOutput< MeshBase >::_serial_only_needed_on_proc_0
protectedinherited

Flag specifying whether this format can be written by only serializing the mesh to processor zero

If this is false (default) the mesh will be serialized to all processors

Definition at line 168 of file mesh_output.h.

◆ line_width

Real libMesh::PostscriptIO::line_width

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 88 of file postscript_io.h.

Referenced by write().

◆ shade_value

Real libMesh::PostscriptIO::shade_value

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 80 of file postscript_io.h.

Referenced by plot_linear_elem().


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