* modify the default outfit of the grids written into a file. See the
* different subclasses and the documentation of the @ref{GridOut}
* class for more details.
- *
- * @author Wolfgang Bangerth, 1998, 2001
*/
namespace GridOutFlags
{
template <>
struct Eps<2> : public EpsFlagsBase
{
+ /**
+ * If this flag is set, then we
+ * place the number of the cell
+ * into the middle of each
+ * cell. The default value is
+ * to not do this.
+ *
+ * The format of the cell
+ * number written is
+ * @p{level.index}.
+ */
+ bool plot_cell_numbers;
+
/**
* Constructor.
*/
const unsigned int size = 300,
const double line_width = 0.5,
const bool color_lines_on_user_flag = false,
- const unsigned int n_boundary_face_points = 2);
+ const unsigned int n_boundary_face_points = 2,
+ const bool plot_cell_numbers = false);
};
/**
const unsigned int size,
const double line_width,
const bool color_lines_on_user_flag,
- const unsigned int n_boundary_face_points)
+ const unsigned int n_boundary_face_points,
+ const bool plot_cell_numbers)
:
EpsFlagsBase(size_type, size, line_width,
color_lines_on_user_flag,
- n_boundary_face_points)
+ n_boundary_face_points),
+ plot_cell_numbers (plot_cell_numbers)
{};
<< "/x {lineto stroke} bind def" << std::endl
<< "/b {0 0 0 setrgbcolor} def" << std::endl
<< "/r {1 0 0 setrgbcolor} def" << std::endl;
+
+ // in 2d, we can also plot cell
+ // numbers, but this requires a
+ // somewhat more lengthy
+ // preamble. please don't ask
+ // me what most of this means,
+ // it is reverse engineered
+ // from what GNUPLOT uses in
+ // its output
+ if ((dim == 2) && (eps_flags_2.plot_cell_numbers == true))
+ {
+ out << ("/R {rmoveto} bind def\n"
+ "/Symbol-Oblique /Symbol findfont [1 0 .167 1 0 0] makefont\n"
+ "dup length dict begin {1 index /FID eq {pop pop} {def} ifelse} forall\n"
+ "currentdict end definefont\n"
+ "/MFshow {{dup dup 0 get findfont exch 1 get scalefont setfont\n"
+ "[ currentpoint ] exch dup 2 get 0 exch rmoveto dup dup 5 get exch 4 get\n"
+ "{show} {stringwidth pop 0 rmoveto}ifelse dup 3 get\n"
+ "{2 get neg 0 exch rmoveto pop} {pop aload pop moveto}ifelse} forall} bind def\n"
+ "/MFwidth {0 exch {dup 3 get{dup dup 0 get findfont exch 1 get scalefont setfont\n"
+ "5 get stringwidth pop add}\n"
+ "{pop} ifelse} forall} bind def\n"
+ "/MCshow { currentpoint stroke m\n"
+ "exch dup MFwidth -2 div 3 -1 roll R MFshow } def\n")
+ << std::endl;
+ };
out << "%%EndProlog" << std::endl
<< std::endl;
<< (line->first - offset) * scale << " m "
<< (line->second - offset) * scale << " x" << std::endl;
+ // finally write the cell numbers
+ // in 2d, if that is desired
+ if ((dim == 2) && (eps_flags_2.plot_cell_numbers == true))
+ {
+ out << "(Helvetica) findfont 140 scalefont setfont"
+ << std::endl;
+
+ typename Triangulation<dim>::active_cell_iterator
+ cell = tria.begin_active (),
+ endc = tria.end ();
+ for (; cell!=endc; ++cell)
+ out << (cell->center()(0)-offset(0))*scale << ' '
+ << (cell->center()(1)-offset(1))*scale
+ << " m" << std::endl
+ << "[ [(Helvetica) 12.0 0.0 true true (" << cell << " )] "
+ << "] -6 MCshow"
+ << std::endl;
+ };
+
out << "showpage" << std::endl;
AssertThrow (out, ExcIO());
<p>
This is the list of changes made after the release of
-<acronym>deal.II</acronym> version 3.3. It is subdivided into changes
+<acronym>deal.II</acronym> version 3.4. It is subdivided into changes
made to the three sub-libraries <a href="#base">base</a>,
<a href="#lac">lac</a>, and <a href="#deal.II">deal.II</a>, as well as
changes to the <a href="#general">general infrastructure,
<h3>deal.II</h3>
<ol>
+ <li> <p>
+ New: For encapsulated postscript of 2d grids, it is now
+ possible to tell the <code class="class">GridOut</code> to
+ write the cell numbers into each cell.
+ <br>
+ (WB 2002/06/04)
+ </p>
</ol>
<hr>
dof_renumbering.exe : dof_renumbering.go $(lib-1d) $(lib-2d) $(lib-3d) $(libraries)
error_estimator.exe : error_estimator.go $(lib-1d) $(lib-2d) $(lib-3d) $(libraries)
gradients.exe : gradients.go $(lib-2d) $(libraries)
+grid_out.exe : grid_out.go $(lib-1d) $(lib-2d) $(lib-3d) $(libraries)
grid_test.exe : grid_test.go $(lib-2d) $(lib-3d) $(libraries)
grid_transform.exe : grid_transform.go $(lib-2d) $(lib-3d) $(libraries)
intergrid_constraints.exe : intergrid_constraints.go $(lib-1d) $(lib-2d) $(lib-3d) $(libraries)
matrices error_estimator intergrid_constraints intergrid_map \
wave-test-3 dof_renumbering support_point_map filtered_matrix \
boundaries sparsity_pattern grid_tools subdomain_ids \
- filtered_iterator
+ filtered_iterator grid_out
############################################################
--- /dev/null
+//---------------------------- constraints.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2002 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//---------------------------- constraints.cc ---------------------------
+
+
+#include <dofs/dof_handler.h>
+#include <grid/tria.h>
+#include <grid/tria_boundary.h>
+#include <grid/tria_boundary_lib.h>
+#include <grid/grid_out.h>
+#include <grid/grid_generator.h>
+#include <base/logstream.h>
+
+#include <fstream>
+
+
+std::ofstream logfile("grid_out.output");
+
+
+template <int dim>
+void test ()
+{
+ Triangulation<dim> tria;
+ if (dim == 2)
+ {
+ static const HyperBallBoundary<dim> x;
+ tria.set_boundary (0, x);
+ GridGenerator::hyper_ball (tria);
+ }
+ else
+ GridGenerator::hyper_cube (tria);
+ tria.refine_global(1);
+
+ GridOut grid_out;
+ GridOutFlags::Eps<2> eps2(GridOutFlags::EpsFlagsBase::width,
+ 300, .5, false, 5, true);
+ grid_out.set_flags (eps2);
+ std::ofstream x("x.eps");
+ grid_out.write_eps (tria, x);
+};
+
+
+int main ()
+{
+ logfile.precision (2);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+
+ test<2> ();
+};
+