<h3>Specific improvements</h3>
<ol>
+ <li> New: There is now a function GridOut::write_vtk() that can
+ write a mesh in VTK format.
+ <br>
+ (Wolfgang Bangerth, 2014/07/14)
+ </li>
+
<li> Fixed: PETSc up to at least version 3.5 has a bug where it does
not zero-initialize the ghost elements of a newly created ghosted
parallel vector. This is now worked around inside deal.II.
#include <deal.II/base/config.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/base/point.h>
+#include <deal.II/base/data_out_base.h>
#include <deal.II/grid/tria.h>
#include <deal.II/fe/mapping.h>
void parse_parameters (ParameterHandler ¶m);
};
+
+ /**
+ * Flags for grid output in Vtk format. These flags are the same as
+ * those declared in DataOutBase::VtkFlags.
+ *
+ * @ingroup output
+ */
+ struct Vtk : public DataOutBase::VtkFlags
+ {};
}
/// write() calls write_svg()
svg,
/// write() calls write_mathgl()
- mathgl
+ mathgl,
+ /// write() calls write_vtk()
+ vtk
};
/**
* colorbar can be drawn to encode the chosen coloring. Moreover, a
* cell label can be added, showing level index, etc.
*
- * @note Yet only implemented for
+ * @note This function is currently only implemented for
* two-dimensional grids in two
* space dimensions.
- *
+ */
+ void write_svg (const Triangulation<2,2> &tria,
+ std::ostream &out) const;
+
+ /**
+ * Declaration of the same function as above for all other dimensions and
+ * space dimensions. This function is not currently implemented and is only
+ * declared to exist to support dimension independent programming.
*/
template <int dim, int spacedim>
void write_svg (const Triangulation<dim,spacedim> &tria,
std::ostream &out) const;
- void write_svg (const Triangulation<2,2> &tria,
- std::ostream &out) const;
/**
* Write triangulation in MathGL script format. To interpret this
*
* To get a handle on the resultant MathGL script within a graphical
* environment an interpreter is needed. A suggestion to start with
- * is <code>mglview</code>, which is bundled with MathGL. With
+ * is <code>mglview</code>, which is bundled with MathGL.
* <code>mglview</code> can interpret and display small-to-medium
* MathGL scripts in a graphical window and enables conversion to
* other formats such as EPS, PNG, JPG, SVG, as well as view/display
* animations. Some minor editing, such as modifying the lighting or
* alpha channels, can also be done.
*
- * @note Not implemented for the codimensional one case.
+ * @note Not implemented for the codimension one case.
*/
template <int dim>
void write_mathgl (const Triangulation<dim> &tria,
std::ostream &out) const;
+ /**
+ * Write triangulation in VTK format.
+ */
+ template <int dim, int spacedim>
+ void write_vtk (const Triangulation<dim,spacedim> &tria,
+ std::ostream &out) const;
+
/**
* Write grid to @p out according to the given data format. This
* function simply calls the appropriate <tt>write_*</tt> function.
void set_flags (const GridOutFlags::MathGL &flags);
/**
- * Provide a function which tells us which
- * suffix with a given output format
- * usually has. At present the following
- * formats are defined:
+ * Set flags for VTK output
+ */
+ void set_flags (const GridOutFlags::Vtk &flags);
+
+ /**
+ * Provide a function that can tell us which
+ * suffix a given output format
+ * usually has. For example, it defines the following mappings:
* <ul>
* <li> @p OpenDX: <tt>.dx</tt>
* <li> @p gnuplot: <tt>.gnuplot</tt>
* <li> @p ucd: <tt>.inp</tt>
* <li> @p eps: <tt>.eps</tt>.
* </ul>
+ * Similar mappings are provided for all implemented formats.
*
* Since this function does not need data from this object, it is
* static and can thus be called without creating an object of this
GridOutFlags::Svg svg_flags;
/**
- * Flags for OpenDX output.
+ * Flags for MathGL output.
*/
GridOutFlags::MathGL mathgl_flags;
+ /**
+ * Flags for VTK output.
+ */
+ GridOutFlags::Vtk vtk_flags;
+
/**
* Write the grid information about faces to @p out. Only those
* faces are printed which are on the boundary and which have a
* 1d. Simply returns zero.
*/
unsigned int n_boundary_faces (const Triangulation<1,1> &tria) const;
+
/**
* Declaration of the specialization of above function for 1d,
* 2sd. Simply returns zero.
*/
unsigned int n_boundary_faces (const Triangulation<1,2> &tria) const;
unsigned int n_boundary_faces (const Triangulation<1,3> &tria) const;
+
/**
* Return the number of lines in the triangulation which have a
* boundary indicator not equal to zero. Only these lines are
:
n_data_sets,
patch->data.n_rows()));
- Assert (patch->data.n_cols() == Utilities::fixed_power<dim>(n_subdivisions+1),
+ Assert ((n_data_sets == 0)
+ ||
+ (patch->data.n_cols() == Utilities::fixed_power<dim>(n_subdivisions+1)),
ExcInvalidDatasetSize (patch->data.n_cols(), n_subdivisions+1));
for (unsigned int i=0; i<patch->data.n_cols(); ++i, ++next_value)
{
draw_bounding_box = param.get_bool ("Draw bounding box");
}
-
} // end namespace GridOutFlags
mathgl_flags = flags;
}
+void GridOut::set_flags (const GridOutFlags::Vtk &flags)
+{
+ vtk_flags = flags;
+}
+
std::string
GridOut::default_suffix (const OutputFormat output_format)
{
return ".svg";
case mathgl:
return ".mathgl";
+ case vtk:
+ return ".vtk";
default:
Assert (false, ExcNotImplemented());
return "";
if (format_name == "mathgl")
return mathgl;
+ if (format_name == "vtk")
+ return vtk;
+
AssertThrow (false, ExcInvalidState ());
// return something weird
return OutputFormat(-1);
std::string GridOut::get_output_format_names ()
{
- return "none|dx|gnuplot|eps|ucd|xfig|msh|svg|mathgl";
+ return "none|dx|gnuplot|eps|ucd|xfig|msh|svg|mathgl|vtk";
}
param.enter_subsection("MathGL");
GridOutFlags::MathGL::declare_parameters(param);
param.leave_subsection();
+
+ param.enter_subsection("Vtk");
+ GridOutFlags::Vtk::declare_parameters(param);
+ param.leave_subsection();
}
param.enter_subsection("MathGL");
mathgl_flags.parse_parameters(param);
param.leave_subsection();
+
+ param.enter_subsection("Vtk");
+ vtk_flags.parse_parameters(param);
+ param.leave_subsection();
}
sizeof(eps_flags_3) +
sizeof(xfig_flags) +
sizeof(svg_flags) +
- sizeof(mathgl_flags));
+ sizeof(mathgl_flags) +
+ sizeof(vtk_flags));
}
}
+
+namespace
+{
+ /**
+ * A function that is able to convert each cell of a triangulation into
+ * a patch that can then be output by the functions in DataOutBase.
+ * This is made particularly simple because the patch only needs to
+ * contain geometry info -- we attach no data at all
+ */
+ template <int dim, int spacedim>
+ std::vector<DataOutBase::Patch<dim,spacedim> >
+ triangulation_to_patches (const Triangulation<dim,spacedim> &triangulation)
+ {
+ std::vector<DataOutBase::Patch<dim,spacedim> > patches;
+ patches.reserve (triangulation.n_active_cells());
+
+ // convert each of the active cells into a patch
+ for (typename Triangulation<dim,spacedim>::active_cell_iterator cell = triangulation.begin_active();
+ cell != triangulation.end(); ++cell)
+ {
+ DataOutBase::Patch<dim,spacedim> patch;
+ patch.n_subdivisions = 1;
+
+ for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
+ patch.vertices[v] = cell->vertex(v);
+
+ // no data
+ patch.data.reinit (0,0);
+
+ patches.push_back (patch);
+ }
+
+ return patches;
+ }
+}
+
+
+
+template <int dim, int spacedim>
+void GridOut::write_vtk (const Triangulation<dim,spacedim> &tria,
+ std::ostream &out) const
+{
+ AssertThrow (out, ExcIO ());
+
+ // convert the cells of the triangulation into a set of patches
+ // and then have them output. since there is no data attached to
+ // the geometry, we also do not have to provide any names, identifying
+ // information, etc.
+ DataOutBase::write_vtk (triangulation_to_patches(tria),
+ std::vector<std::string>(),
+ std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> >(),
+ vtk_flags,
+ out);
+
+ AssertThrow (out, ExcIO ());
+}
+
+
+
unsigned int GridOut::n_boundary_faces (const Triangulation<1> &) const
{
return 0;
case mathgl:
write_mathgl (tria, out);
return;
+
+ case vtk:
+ write_vtk (tria, out);
+ return;
}
Assert (false, ExcInternalError());
template <int dim, int spacedim>
-void GridOut::write (
- const Triangulation<dim, spacedim> &tria,
- std::ostream &out,
- const Mapping<dim,spacedim> *mapping) const
+void GridOut::write (const Triangulation<dim, spacedim> &tria,
+ std::ostream &out,
+ const Mapping<dim,spacedim> *mapping) const
{
write(tria, out, default_format, mapping);
}
(const Triangulation<deal_II_dimension> &,
std::ostream &,
const Mapping<deal_II_dimension,deal_II_dimension> *) const;
+ template void GridOut::write_vtk
+ (const Triangulation<deal_II_dimension>&,
+ std::ostream&) const;
+
template void GridOut::write<deal_II_dimension>
(const Triangulation<deal_II_dimension> &,
std::ostream &, const OutputFormat,
(const Triangulation<deal_II_dimension,deal_II_dimension+1>&,
std::ostream&,
const Mapping<deal_II_dimension,deal_II_dimension+1>*) const;
+ template void GridOut::write_vtk
+ (const Triangulation<deal_II_dimension,deal_II_dimension+1>&,
+ std::ostream&) const;
#endif
#if deal_II_dimension == 3
template void GridOut::write_ucd