// $Id$
// Version: $Name$
//
-// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008 by the deal.II authors
+// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2010 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
* Write cells.
*/
bool write_cells;
-
+
/**
* Write faces.
*/
bool write_faces;
-
+
/**
* Write field with diameters.
*/
bool write_diameter;
-
+
/**
* Write field with area/volume.
*/
- bool write_measure;
-
+ bool write_measure;
+
/**
* Write all faces, including
* interior faces. If
* boundary faces are written.
*/
bool write_all_faces;
-
+
/**
* Constructor.
*/
* ParameterHandler.
*/
void parse_parameters (ParameterHandler& param);
- };
-
+ };
+
/**
* Flags describing the details of
* output in MSH format.
*
* @ingroup output
*/
- struct Msh
+ struct Msh
{
/**
* When writing a mesh, write boundary
* It is not necessary if you only want
* to write the triangulation to
* view or print it.
- *
+ *
* This is used only if
* <tt>dim==3</tt>, and ignored
* in all other cases.
* Default: @p false.
*/
bool write_lines;
-
+
/**
* Constructor.
*/
- Msh (const bool write_faces = false,
+ Msh (const bool write_faces = false,
const bool write_lines = false);
/**
* Declare parameters in
*/
void parse_parameters (ParameterHandler& param);
};
-
-
+
+
/**
* Flags describing the details of
* output in UCD format.
*
* @ingroup output
*/
- struct Ucd
+ struct Ucd
{
/**
* Write a comment at the
* Default: <code>false</code>.
*/
bool write_preamble;
-
+
/**
* When writing a mesh, write boundary
* faces explicitly if their boundary
* Default: @p false.
*/
bool write_faces;
-
+
/**
* When writing a mesh, write boundary
* lines explicitly if their boundary
* Default: @p false.
*/
bool write_lines;
-
+
/**
* Constructor.
*/
*/
void parse_parameters (ParameterHandler& param);
};
-
-
+
+
/**
* Flags describing the details of
* output in GNUPLOT format.
* <code>n_boundary_face_points>0</code>.
*/
bool curved_inner_cells;
-
+
/**
* Constructor.
*/
enum SizeType {
width, height
};
-
+
/**
* See above. Default is @p width.
*/
SizeType size_type;
-
+
/**
* Width or height of the output
* as given in postscript units
* Default is 300.
*/
unsigned int size;
-
+
/**
* Width of a line in postscript
* units. Default is 0.5.
*/
double line_width;
-
+
/**
* Should lines with a set
* @p user_flag be drawn in a
* different color (red)?
*/
bool color_lines_on_user_flag;
-
+
/**
* This is the number of
* points on a boundary face,
* mapping.
*/
unsigned int n_boundary_face_points;
-
+
/**
* Should lines be colored
* according to their
* from blue to red.
*/
bool color_lines_level;
-
+
/**
* Constructor.
*/
*/
void parse_parameters (ParameterHandler& param);
};
-
-
+
+
/**
* Flags describing the details of
* output for encapsulated postscript
template <int dim>
struct Eps
{};
-
+
/**
* Flags specific to the output of
* grids in one space dimensions.
* @ingroup output
*/
template <>
- struct Eps<1> : public EpsFlagsBase
+ struct Eps<1> : public EpsFlagsBase
{
/**
* Constructor.
void parse_parameters (ParameterHandler& param);
};
-
+
/**
* Flags specific to the output of
* grids in two space dimensions.
* @ingroup output
*/
template <>
- struct Eps<2> : public EpsFlagsBase
+ struct Eps<2> : public EpsFlagsBase
{
/**
* If this flag is set, then we
* flag. Default is @p false.
*/
bool write_vertex_numbers;
-
+
/**
* Constructor.
*/
*/
void parse_parameters (ParameterHandler& param);
};
-
+
/**
* Flags specific to the output of
* grids in three space dimensions.
* of 60.
*/
double azimut_angle;
-
+
/**
* Angle by which the viewers
* position projected onto the
* of 30.
*/
double turn_angle;
-
+
/**
* Constructor.
*/
* is forwarded to XFig.
*/
int line_style;
-
+
/**
* Thickness of border lines of
* polygons. Default is 1.
* meshes.
*/
int line_thickness;
-
+
/**
* Style for drawing lines at the
* boundary. Defaults to solid (0).
*/
int boundary_style;
-
+
/**
* Thickness of boundary
* lines. Default is 3.
*/
int boundary_thickness;
-
+
/**
* Constructor.
*/
*/
void parse_parameters (ParameterHandler& param);
};
-
+
}
* @ingroup output
* @author Wolfgang Bangerth, Guido Kanschat, Luca Heltai, 1999, 2003, 2006; postscript format based on an implementation by Stefan Nauber, 1999
*/
-class GridOut
+class GridOut
{
public:
/**
/// write() calls write_ucd()
ucd,
/// write() calls write_xfig()
- xfig,
+ xfig,
/// write() calls write_msh()
msh
};
-
+
+ /**
+ * Constructor.
+ */
+ GridOut ();
+
/**
* Write triangulation in OpenDX
* format.
* together with their level and
* their material id or boundary
* indicator, resp.
- *
+ *
* Not implemented for the
* codimension one case.
*/
void write_gnuplot (const Triangulation<dim> &tria,
std::ostream &out,
const Mapping<dim> *mapping=0) const;
-
+
/**
* Write the triangulation in the
* msh format.
*
- * Msh
+ * Msh
* is the format used by Gmsh and it is
* described in the gmsh
* user's guide. Besides the
* the extension of the picture,
* of which the default is 300.
*
- * The flag
+ * The flag
* @p color_lines_on_user_flag
* allows to draw lines with the
* @p user_flag set to be drawn in
void write_xfig (const Triangulation<dim> &tria,
std::ostream &out,
const Mapping<dim> *mapping=0) const;
-
+
/**
* Write grid to @p out according
* to the given data format. This
std::ostream &out,
const OutputFormat output_format,
const Mapping<dim,spacedim> *mapping=0) const;
-
+
/**
* Write mesh in default format
* set by ParameterHandler.
void write (const Triangulation<dim,spacedim> &tria,
std::ostream &out,
const Mapping<dim,spacedim> *mapping=0) const;
-
+
/**
* Set flags for DX output
*/
* ParameterHandler.
*/
std::string default_suffix () const;
-
+
/**
* Return the @p OutputFormat value
* corresponding to the given string. If
* ParameterHandler.
*/
static void declare_parameters (ParameterHandler& param);
-
+
/**
* Parse parameters of
* ParameterHandler.
*/
void parse_parameters (ParameterHandler& param);
-
+
/**
* Determine an estimate for the
* memory consumption (in bytes)
* by a ParameterHandler.
*/
OutputFormat default_format;
-
+
/**
* Flags for OpenDX output.
*/
GridOutFlags::DX dx_flags;
-
+
/**
* Flags for GMSH output. Can be
* changed by using the
* @p set_flags function.
*/
GridOutFlags::Msh msh_flags;
-
+
/**
* Flags for UCD output. Can be
* changed by using the
* function.
*/
GridOutFlags::Eps<3> eps_flags_3;
-
+
/**
* Flags used for XFig output.
*/
GridOutFlags::XFig xfig_flags;
-
+
/**
* Write the grid information about
* faces to @p out. Only those faces
void write_msh_faces (const Triangulation<1,2> &tria,
const unsigned int starting_index,
std::ostream &out) const;
-
+
/**
* Write the grid information about
* lines to @p out. Only those lines
void write_msh_lines (const Triangulation<1,1> &tria,
const unsigned int starting_index,
std::ostream &out) const;
-
+
/**
* Declaration of the specialization
* of above function for 1d, 2sd. Does
* the function for <tt>dim==1</tt>. Bad luck.
*/
-
+
template <int dim, int spacedim>
void write_ucd_faces (const Triangulation<dim,spacedim> &tria,
const unsigned int starting_index,
* the function for <tt>dim==1/2</tt>. Bad luck.
*/
-
+
template <int dim, int spacedim>
void write_ucd_lines (const Triangulation<dim,spacedim> &tria,
const unsigned int starting_index,
void write_ucd_lines (const Triangulation<2,3> &tria,
const unsigned int starting_index,
std::ostream &out) const;
-
+
/**
* Return the number of faces in the
* triangulation which have a boundary
// $Id$
// Version: $Name$
//
-// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008 by the deal.II authors
+// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2010 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
DEAL_II_NAMESPACE_OPEN
+GridOut::GridOut ()
+ :
+ default_format (none)
+{}
+
+
#if deal_II_dimension == 1
void GridOut::write_dx (const Triangulation<dim> &tria,
std::ostream &out) const
{
-//TODO:[GK] allow for boundary faces only
+//TODO:[GK] allow for boundary faces only
Assert(dx_flags.write_all_faces, ExcNotImplemented());
AssertThrow (out, ExcIO());
// Copied and adapted from write_ucd
const std::vector<Point<dim> > &vertices = tria.get_vertices();
const std::vector<bool> &vertex_used = tria.get_used_vertices();
-
+
const unsigned int n_vertices = tria.n_used_vertices();
// vertices are implicitly numbered from 0 to
typename Triangulation<dim>::active_cell_iterator cell;
const typename Triangulation<dim>::active_cell_iterator endc=tria.end();
-
+
// write the vertices
out << "object \"vertices\" class array type float rank 1 shape " << dim
<< " items " << n_vertices << " data follows"
<< '\n';
-
+
for (unsigned int i=0; i<vertices.size(); ++i)
if (vertex_used[i])
{
out << '\t' << vertices[i] << '\n';
};
-
+
// write cells or faces
const bool write_cells = dx_flags.write_cells;
const bool write_faces = (dim>1) ? dx_flags.write_faces : false;
-
+
const unsigned int n_cells = tria.n_active_cells();
const unsigned int n_faces = tria.n_active_cells()
* GeometryInfo<dim>::faces_per_cell;
const unsigned int n_vertices_per_cell = GeometryInfo<dim>::vertices_per_cell;
const unsigned int n_vertices_per_face = GeometryInfo<dim>::vertices_per_face;
-
+
if (write_cells)
{
out << "object \"cells\" class array type int rank 1 shape "
<< n_vertices_per_cell
<< " items " << n_cells << " data follows" << '\n';
-
+
for (cell = tria.begin_active(); cell != endc; ++cell)
{
for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
<< "attribute \"ref\" string \"positions\"" << '\n' << '\n';
// Additional cell information
-
+
out << "object \"material\" class array type int rank 0 items "
<< n_cells << " data follows" << '\n';
for (cell = tria.begin_active(); cell != endc; ++cell)
out << ' ' << (unsigned int)cell->material_id();
out << '\n'
<< "attribute \"dep\" string \"connections\"" << '\n' << '\n';
-
+
out << "object \"level\" class array type int rank 0 items "
<< n_cells << " data follows" << '\n';
for (cell = tria.begin_active(); cell != endc; ++cell)
out << ' ' << cell->level();
out << '\n'
<< "attribute \"dep\" string \"connections\"" << '\n' << '\n';
-
+
if (dx_flags.write_measure)
{
out << "object \"measure\" class array type float rank 0 items "
out << '\n'
<< "attribute \"dep\" string \"connections\"" << '\n' << '\n';
}
-
+
if (dx_flags.write_diameter)
{
out << "object \"diameter\" class array type float rank 0 items "
<< "attribute \"dep\" string \"connections\"" << '\n' << '\n';
}
}
-
+
if (write_faces)
{
out << "object \"faces\" class array type int rank 1 shape "
for (unsigned int f=0;f<GeometryInfo<dim>::faces_per_cell;++f)
{
typename Triangulation<dim>::face_iterator face = cell->face(f);
-
+
for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
out << '\t' << renumber[face->vertex_index(GeometryInfo<dim-1>::dx_to_deal[v])];
out << '\n';
if (dim==3) out << "quads";
out << "\"" << '\n'
<< "attribute \"ref\" string \"positions\"" << '\n' << '\n';
-
+
// Additional face information
-
+
out << "object \"boundary\" class array type int rank 0 items "
<< n_faces << " data follows" << '\n';
for (cell = tria.begin_active(); cell != endc; ++cell)
out << '\n';
}
out << "attribute \"dep\" string \"connections\"" << '\n' << '\n';
-
+
if (dx_flags.write_measure)
{
out << "object \"face measure\" class array type float rank 0 items "
}
}
-
+
// Write additional face information
-
+
if (write_faces)
{
-
+
}
else
{
if (dx_flags.write_diameter)
out << "component \"diameter\" value \"face diameter\"" << '\n';
}
-
+
out << '\n'
<< "object \"grid data\" class group" << '\n';
if (write_cells)
out << "member \"cells\" value \"cell data\"" << '\n';
if (write_faces)
out << "member \"faces\" value \"face data\"" << '\n';
- out << "end" << '\n';
+ out << "end" << '\n';
// make sure everything now gets to
// disk
// used.
const std::vector<Point<spacedim> > &vertices = tria.get_vertices();
const std::vector<bool> &vertex_used = tria.get_used_vertices();
-
+
const unsigned int n_vertices = tria.n_used_vertices();
typename Triangulation<dim,spacedim>::active_cell_iterator cell=tria.begin_active();
// Write Header
// The file format is:
- /*
+ /*
+
-
$NOD
number-of-nodes
node-number x-coord y-coord z-coord
...
$ENDELM
*/
- out << "$NOD" << std::endl
+ out << "$NOD" << std::endl
<< n_vertices << std::endl;
// actually write the vertices.
out << "$ENDNOD" << std::endl
<< "$ELM" << std::endl
<< tria.n_active_cells() + ( (msh_flags.write_faces ?
- n_boundary_faces(tria) :0 ) +
- ( msh_flags.write_lines ?
+ n_boundary_faces(tria) :0 ) +
+ ( msh_flags.write_lines ?
n_boundary_lines(tria) : 0 ) ) << std::endl;
/*
{
out << cell_index << ' ' << elm_type << ' '
<< static_cast<unsigned int>(cell->material_id()) << ' '
- << cell->subdomain_id() << ' '
+ << cell->subdomain_id() << ' '
<< GeometryInfo<dim>::vertices_per_cell << ' ';
// Vertex numbering follows UCD conventions.
// write faces and lines with
// non-zero boundary indicator
- if (msh_flags.write_faces)
+ if (msh_flags.write_faces)
write_msh_faces (tria, cell_index, out);
if (msh_flags.write_lines)
write_msh_lines (tria, cell_index, out);
// variables destroyed after
// use
std::time_t time1= std::time (0);
- std::tm *time = std::localtime(&time1);
+ std::tm *time = std::localtime(&time1);
out << "# This file was generated by the deal.II library." << '\n'
<< "# Date = "
<< time->tm_year+1900 << "/"
out << n_vertices << ' '
<< tria.n_active_cells() + ( (ucd_flags.write_faces ?
n_boundary_faces(tria) : 0) +
- (ucd_flags.write_lines ?
+ (ucd_flags.write_lines ?
n_boundary_lines(tria) : 0) )
<< " 0 0 0" // no data
<< '\n';
out << cell_index << ' '
<< static_cast<unsigned int>(cell->material_id())
<< ' ';
- switch (dim)
+ switch (dim)
{
case 1: out << "line "; break;
case 2: out << "quad "; break;
// AVS Developer's Guide, Release 4,
// May, 1992, p. E6
//
- // note: vertex numbers are 1-base
+ // note: vertex numbers are 1-base
for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
++vertex)
out << cell->vertex_index(GeometryInfo<dim>::ucd_to_deal[vertex])+1 << ' ';
// write faces and lines with
// non-zero boundary indicator
- if (ucd_flags.write_faces)
+ if (ucd_flags.write_faces)
write_ucd_faces (tria, cell_index, out);
if (ucd_flags.write_lines)
write_ucd_lines (tria, cell_index, out);
out << '\t' << val;
}
out << '\n';
- }
+ }
}
}
}
for (face=tria.begin_active_face(), endf=tria.end_face();
face != endf; ++face)
if (face->at_boundary() &&
- (face->boundary_indicator() != 0))
+ (face->boundary_indicator() != 0))
{
out << index << ' ';
- switch (dim)
+ switch (dim)
{
case 2: out << 1 << ' '; break;
case 3: out << 3 << ' '; break;
default:
Assert (false, ExcNotImplemented());
};
- out << static_cast<unsigned int>(face->boundary_indicator())
- << ' '
- << static_cast<unsigned int>(face->boundary_indicator())
+ out << static_cast<unsigned int>(face->boundary_indicator())
+ << ' '
+ << static_cast<unsigned int>(face->boundary_indicator())
<< ' ' << GeometryInfo<dim>::vertices_per_face;
// note: vertex numbers are 1-base
for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
- out << ' '
+ out << ' '
<< face->vertex_index(GeometryInfo<dim-1>::ucd_to_deal[vertex])+1;
out << '\n';
++index;
- };
+ };
}
for (line=tria.begin_active_line(), endl=tria.end_line();
line != endl; ++line)
if (line->at_boundary() &&
- (line->boundary_indicator() != 0))
+ (line->boundary_indicator() != 0))
{
out << index << " 1 ";
- out << static_cast<unsigned int>(line->boundary_indicator())
- << ' '
- << static_cast<unsigned int>(line->boundary_indicator())
+ out << static_cast<unsigned int>(line->boundary_indicator())
+ << ' '
+ << static_cast<unsigned int>(line->boundary_indicator())
<< " 2 ";
// note: vertex numbers are 1-base
for (unsigned int vertex=0; vertex<2; ++vertex)
- out << ' '
+ out << ' '
<< line->vertex_index(GeometryInfo<dim-2>::ucd_to_deal[vertex])+1;
out << '\n';
++index;
- };
+ };
}
for (face=tria.begin_active_face(), endf=tria.end_face();
face != endf; ++face)
if (face->at_boundary() &&
- (face->boundary_indicator() != 0))
+ (face->boundary_indicator() != 0))
{
out << index << " "
<< static_cast<unsigned int>(face->boundary_indicator())
<< " ";
- switch (dim)
+ switch (dim)
{
case 2: out << "line "; break;
case 3: out << "quad "; break;
out << '\n';
++index;
- };
+ };
}
for (line=tria.begin_active_line(), endl=tria.end_line();
line != endl; ++line)
if (line->at_boundary() &&
- (line->boundary_indicator() != 0))
+ (line->boundary_indicator() != 0))
{
out << index << " "
<< static_cast<unsigned int>(line->boundary_indicator())
out << '\n';
++index;
- };
+ };
}
{
boundary_points.resize(n_points);
boundary_points[0][0]=0;
- boundary_points[n_points-1][0]=1;
+ boundary_points[n_points-1][0]=1;
for (unsigned int i=1; i<n_points-1; ++i)
boundary_points[i](0)= 1.*i/(n_points-1);
}
if (q_projector != 0)
- delete q_projector;
+ delete q_projector;
// make sure everything now gets to
// disk
{
boundary_points.resize(n_points);
boundary_points[0][0]=0;
- boundary_points[n_points-1][0]=1;
+ boundary_points[n_points-1][0]=1;
for (unsigned int i=1; i<n_points-1; ++i)
boundary_points[i](0)= 1.*i/(n_points-1);
(dim==3 ?
static_cast<const GridOutFlags::EpsFlagsBase&>(eps_flags_3) :
*static_cast<const GridOutFlags::EpsFlagsBase*>(0)));
-
+
AssertThrow (out, ExcIO());
const unsigned int n_points = eps_flags_base.n_boundary_face_points;
-
+
// make up a list of lines by which
// we will construct the triangulation
//
// actual output dimension independent
// again
LineList line_list;
-
+
switch (dim)
{
case 1:
Assert(false, ExcInternalError());
break;
};
-
+
case 2:
{
typename Triangulation<dim>::active_cell_iterator
// treat all other
// lines
if (!line->has_children() &&
- (mapping==0 || !line->at_boundary()))
+ (mapping==0 || !line->at_boundary()))
// one would expect
// make_pair(line->vertex(0),
- // line->vertex(1))
+ // line->vertex(1))
// here, but that is
// not dimension
// independent, since
line->user_flag_set(),
cell->level()));
}
-
+
// next if we are to treat
// curved boundaries
// specially, then add lines
// project them onto the
// faces of a unit cell
std::vector<Point<dim-1> > boundary_points (n_points);
-
+
for (unsigned int i=0; i<n_points; ++i)
boundary_points[i](0) = 1.*(i+1)/(n_points+1);
-
+
Quadrature<dim-1> quadrature (boundary_points);
Quadrature<dim> q_projector (QProjector<dim>::project_to_all_faces(quadrature));
-
+
// next loop over all
// boundary faces and
// generate the info from
// projection on the plane perpendicular
// to the direction of sight
- // direction of view equals the unit
+ // direction of view equals the unit
// vector of the position of the
// spectator to the origin.
//
line->user_flag_set(),
cell->level()));
}
-
+
break;
}
y_max = std::max (y_max, line->first(1));
y_max = std::max (y_max, line->second(1));
- max_level = std::max (max_level, line->level);
+ max_level = std::max (max_level, line->level);
};
// scale in x-direction such that
// now write preamble
- if (true)
+ if (true)
{
// block this to have local
// variables destroyed after
// use
std::time_t time1= std::time (0);
- std::tm *time = std::localtime(&time1);
+ std::tm *time = std::localtime(&time1);
out << "%!PS-Adobe-2.0 EPSF-1.2" << '\n'
<< "%%Title: deal.II Output" << '\n'
<< "%%Creator: the deal.II library" << '\n'
- << "%%Creation Date: "
+ << "%%Creation Date: "
<< time->tm_year+1900 << "/"
<< time->tm_mon+1 << "/"
<< time->tm_mday << " - "
{
case none:
return;
-
+
case dx:
write_dx (tria, out);
return;
-
+
case ucd:
write_ucd (tria, out);
return;
-
+
case gnuplot:
write_gnuplot (tria, out, mapping);
return;
-
+
case eps:
write_eps (tria, out, mapping);
return;
-
+
case xfig:
write_xfig (tria, out, mapping);
return;
-
+
case msh:
write_msh (tria, out);
return;
}
-
+
Assert (false, ExcInternalError());
}