From 4ec93114dfbf6a9a564b043603847243fe80448c Mon Sep 17 00:00:00 2001 From: bangerth Date: Mon, 14 Jul 2014 22:20:35 +0000 Subject: [PATCH] Provide a function GridOut::write_vtk(). git-svn-id: https://svn.dealii.org/trunk@33172 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/doc/news/changes.h | 6 ++ deal.II/include/deal.II/grid/grid_out.h | 60 ++++++++++++---- deal.II/source/base/data_out_base.cc | 4 +- deal.II/source/grid/grid_out.cc | 94 +++++++++++++++++++++++-- deal.II/source/grid/grid_out.inst.in | 7 ++ 5 files changed, 151 insertions(+), 20 deletions(-) diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 265a0f3263..ab287d9439 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -172,6 +172,12 @@ inconvenience this causes.

Specific improvements

    +
  1. New: There is now a function GridOut::write_vtk() that can + write a mesh in VTK format. +
    + (Wolfgang Bangerth, 2014/07/14) +
  2. +
  3. 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. diff --git a/deal.II/include/deal.II/grid/grid_out.h b/deal.II/include/deal.II/grid/grid_out.h index a72565cfd7..29e7742eeb 100644 --- a/deal.II/include/deal.II/grid/grid_out.h +++ b/deal.II/include/deal.II/grid/grid_out.h @@ -22,6 +22,7 @@ #include #include #include +#include #include #include @@ -743,6 +744,15 @@ namespace GridOutFlags 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 + {}; } @@ -851,7 +861,9 @@ public: /// write() calls write_svg() svg, /// write() calls write_mathgl() - mathgl + mathgl, + /// write() calls write_vtk() + vtk }; /** @@ -1066,17 +1078,22 @@ public: * 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 void write_svg (const Triangulation &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 @@ -1084,19 +1101,26 @@ public: * * To get a handle on the resultant MathGL script within a graphical * environment an interpreter is needed. A suggestion to start with - * is mglview, which is bundled with MathGL. With + * is mglview, which is bundled with MathGL. * mglview 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 void write_mathgl (const Triangulation &tria, std::ostream &out) const; + /** + * Write triangulation in VTK format. + */ + template + void write_vtk (const Triangulation &tria, + std::ostream &out) const; + /** * Write grid to @p out according to the given data format. This * function simply calls the appropriate write_* function. @@ -1168,16 +1192,21 @@ public: 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: *
      *
    • @p OpenDX: .dx *
    • @p gnuplot: .gnuplot *
    • @p ucd: .inp *
    • @p eps: .eps. *
    + * 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 @@ -1300,10 +1329,15 @@ private: 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 @@ -1544,12 +1578,14 @@ private: * 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 diff --git a/deal.II/source/base/data_out_base.cc b/deal.II/source/base/data_out_base.cc index 39a23d8558..b0ee12eb5c 100644 --- a/deal.II/source/base/data_out_base.cc +++ b/deal.II/source/base/data_out_base.cc @@ -407,7 +407,9 @@ namespace DataOutBase : n_data_sets, patch->data.n_rows())); - Assert (patch->data.n_cols() == Utilities::fixed_power(n_subdivisions+1), + Assert ((n_data_sets == 0) + || + (patch->data.n_cols() == Utilities::fixed_power(n_subdivisions+1)), ExcInvalidDatasetSize (patch->data.n_cols(), n_subdivisions+1)); for (unsigned int i=0; idata.n_cols(); ++i, ++next_value) diff --git a/deal.II/source/grid/grid_out.cc b/deal.II/source/grid/grid_out.cc index d4415d620a..abbf96cebe 100644 --- a/deal.II/source/grid/grid_out.cc +++ b/deal.II/source/grid/grid_out.cc @@ -392,7 +392,6 @@ namespace GridOutFlags { draw_bounding_box = param.get_bool ("Draw bounding box"); } - } // end namespace GridOutFlags @@ -468,6 +467,11 @@ void GridOut::set_flags (const GridOutFlags::MathGL &flags) mathgl_flags = flags; } +void GridOut::set_flags (const GridOutFlags::Vtk &flags) +{ + vtk_flags = flags; +} + std::string GridOut::default_suffix (const OutputFormat output_format) { @@ -491,6 +495,8 @@ GridOut::default_suffix (const OutputFormat output_format) return ".svg"; case mathgl: return ".mathgl"; + case vtk: + return ".vtk"; default: Assert (false, ExcNotImplemented()); return ""; @@ -537,6 +543,9 @@ GridOut::parse_output_format (const std::string &format_name) if (format_name == "mathgl") return mathgl; + if (format_name == "vtk") + return vtk; + AssertThrow (false, ExcInvalidState ()); // return something weird return OutputFormat(-1); @@ -546,7 +555,7 @@ GridOut::parse_output_format (const std::string &format_name) 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"; } @@ -586,6 +595,10 @@ GridOut::declare_parameters(ParameterHandler ¶m) param.enter_subsection("MathGL"); GridOutFlags::MathGL::declare_parameters(param); param.leave_subsection(); + + param.enter_subsection("Vtk"); + GridOutFlags::Vtk::declare_parameters(param); + param.leave_subsection(); } @@ -624,6 +637,10 @@ GridOut::parse_parameters(ParameterHandler ¶m) param.enter_subsection("MathGL"); mathgl_flags.parse_parameters(param); param.leave_subsection(); + + param.enter_subsection("Vtk"); + vtk_flags.parse_parameters(param); + param.leave_subsection(); } @@ -640,7 +657,8 @@ GridOut::memory_consumption () const sizeof(eps_flags_3) + sizeof(xfig_flags) + sizeof(svg_flags) + - sizeof(mathgl_flags)); + sizeof(mathgl_flags) + + sizeof(vtk_flags)); } @@ -2367,6 +2385,65 @@ void GridOut::write_mathgl (const Triangulation &tria, } + +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 + std::vector > + triangulation_to_patches (const Triangulation &triangulation) + { + std::vector > patches; + patches.reserve (triangulation.n_active_cells()); + + // convert each of the active cells into a patch + for (typename Triangulation::active_cell_iterator cell = triangulation.begin_active(); + cell != triangulation.end(); ++cell) + { + DataOutBase::Patch patch; + patch.n_subdivisions = 1; + + for (unsigned int v=0; v::vertices_per_cell; ++v) + patch.vertices[v] = cell->vertex(v); + + // no data + patch.data.reinit (0,0); + + patches.push_back (patch); + } + + return patches; + } +} + + + +template +void GridOut::write_vtk (const Triangulation &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::vector >(), + vtk_flags, + out); + + AssertThrow (out, ExcIO ()); +} + + + unsigned int GridOut::n_boundary_faces (const Triangulation<1> &) const { return 0; @@ -3736,6 +3813,10 @@ void GridOut::write (const Triangulation &tria, case mathgl: write_mathgl (tria, out); return; + + case vtk: + write_vtk (tria, out); + return; } Assert (false, ExcInternalError()); @@ -3743,10 +3824,9 @@ void GridOut::write (const Triangulation &tria, template -void GridOut::write ( - const Triangulation &tria, - std::ostream &out, - const Mapping *mapping) const +void GridOut::write (const Triangulation &tria, + std::ostream &out, + const Mapping *mapping) const { write(tria, out, default_format, mapping); } diff --git a/deal.II/source/grid/grid_out.inst.in b/deal.II/source/grid/grid_out.inst.in index 1da4845a0b..0b84908ff4 100644 --- a/deal.II/source/grid/grid_out.inst.in +++ b/deal.II/source/grid/grid_out.inst.in @@ -49,6 +49,10 @@ for (deal_II_dimension : DIMENSIONS) (const Triangulation &, std::ostream &, const Mapping *) const; + template void GridOut::write_vtk + (const Triangulation&, + std::ostream&) const; + template void GridOut::write (const Triangulation &, std::ostream &, const OutputFormat, @@ -69,6 +73,9 @@ for (deal_II_dimension : DIMENSIONS) (const Triangulation&, std::ostream&, const Mapping*) const; + template void GridOut::write_vtk + (const Triangulation&, + std::ostream&) const; #endif #if deal_II_dimension == 3 template void GridOut::write_ucd -- 2.39.5