From: bangerth
Date: Sat, 1 Feb 2014 13:18:44 +0000 (+0000)
Subject: Remove the awkward private inheritance from DataOutBase to DataOutInterface that...
X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=62af64704c05f22a436f3bbdad767965436741bb;p=dealii-svn.git
Remove the awkward private inheritance from DataOutBase to DataOutInterface that has tripped up virtually every compiler we have ever tried. Instead, make DataOutBase a namespace.
git-svn-id: https://svn.dealii.org/trunk@32367 0785d39b-7218-0410-832d-ea1e28bc413d
---
diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h
index 73d5cb20e9..384bfde3a7 100644
--- a/deal.II/doc/news/changes.h
+++ b/deal.II/doc/news/changes.h
@@ -40,6 +40,28 @@ inconvenience this causes.
+ - Changed: The various classes generating graphical output, such
+ as DataOut or DataOutStack, are all derived from a common interface
+ class DataOutInterface which, in turn was derived from DataOutBase
+ through private inheritance. Because we frequently also
+ access the (public) members of this private base class this has tripped
+ up most every compiler we know of at one point or another. Furthermore,
+ because DataOutBase was a class that only defined static member functions
+ and had not member variables, there was really no reason for this
+ construct.
+
+ For these reasons, DataOutBase is now just a regular namespace and the
+ inheritance is gone. For the most part, this should not lead to any
+ incompatibilities except in cases where you accessed members of
+ DataOutBase through their derived classes. For example, it was possible
+ to write DataOut@<2@>::Patch@<2,2@>
even though the
+ Patch
class is actually declared in DataOutBase. Since
+ the inheritance is now gone, this is no longer possible and one
+ actually has to write DataOutBase::Patch instead. Using this form
+ turns out to be compatible also with older versions of deal.II.
+
+ (Wolfgang Bangerth, 2014/02/01)
+
diff --git a/deal.II/include/deal.II/base/data_out_base.h b/deal.II/include/deal.II/base/data_out_base.h
index 0d9ee2bbf3..ac5dc7dffe 100644
--- a/deal.II/include/deal.II/base/data_out_base.h
+++ b/deal.II/include/deal.II/base/data_out_base.h
@@ -44,7 +44,6 @@ DEAL_II_NAMESPACE_OPEN
class ParameterHandler;
class XDMFEntry;
-class DataOutFilter;
/**
* This is a base class for output of data on meshes of very general
@@ -216,9 +215,8 @@ class DataOutFilter;
* @ingroup output
* @author Wolfgang Bangerth, Guido Kanschat 1999, 2000, 2001, 2002, 2005, 2006.
*/
-class DataOutBase
+namespace DataOutBase
{
-public:
/**
* Data structure describing a patch of data in dim space
* dimensions.
@@ -2020,151 +2018,13 @@ public:
<< " for output");
//@}
-private:
- /**
- * Write the coordinates of nodes in the desired format.
- */
- template
- static void write_nodes (const std::vector > &patches,
- STREAM &out);
-
- /**
- * Write the node numbers of a cell in the desired format.
- */
- template
- static void write_cells (const std::vector > &patches,
- STREAM &out);
-
- /**
- * Write data in the desired format.
- */
- template
- static void write_data (const std::vector > &patches,
- const unsigned int n_data_sets,
- const bool double_precision,
- STREAM &out);
-
-
- /**
- * This function projects a three-dimensional point (Point<3> point)
- * onto a two-dimensional image plane, specified by the position of
- * the camera viewing system (Point<3> camera_position), camera
- * direction (Point<3> camera_position), camera horizontal (Point<3>
- * camera_horizontal, necessary for the correct alignment of the
- * later images), and the focus of the camera (float camera_focus).
- *
- * For SVG output.
- */
- static Point<2> svg_project_point(Point<3> point,
- Point<3> camera_position,
- Point<3> camera_direction,
- Point<3> camera_horizontal,
- float camera_focus);
- /**
- * Function to compute the gradient parameters for a triangle with
- * given values for the vertices.
- *
- * Used for svg output.
- */
- static Point<6> svg_get_gradient_parameters(Point<3> points[]);
-
- /**
- * Class holding the data of one cell of a patch in two space
- * dimensions for output. It is the projection of a cell in
- * three-dimensional space (two coordinates, one height value) to
- * the direction of sight.
- */
- class SvgCell
- {
- public:
-
- // Center of the cell (three-dimensional)
- Point<3> center;
-
- /**
- * Vector of vertices of this cell (three-dimensional)
- */
- Point<3> vertices[4];
-
- /**
- * Depth into the picture, which is defined as the distance from
- * an observer at an the origin in direction of the line of sight.
- */
- float depth;
-
- /**
- * Vector of vertices of this cell (projected, two-dimensional).
- */
- Point<2> projected_vertices[4];
-
- // Center of the cell (projected, two-dimensional)
- Point<2> projected_center;
-
- /**
- * Comparison operator for sorting.
- */
- bool operator < (const SvgCell &) const;
- };
-
-
- /**
- * Class holding the data of one cell of a patch in two space
- * dimensions for output. It is the projection of a cell in
- * three-dimensional space (two coordinates, one height value) to
- * the direction of sight.
- */
- class EpsCell2d
- {
- public:
-
- /**
- * Vector of vertices of this cell.
- */
- Point<2> vertices[4];
-
- /**
- * Data value from which the actual colors will be computed by the
- * colorization function stated in the EpsFlags class.
- */
- float color_value;
-
- /**
- * Depth into the picture, which is defined as the distance from
- * an observer at an the origin in direction of the line of sight.
- */
- float depth;
-
- /**
- * Comparison operator for sorting.
- */
- bool operator < (const EpsCell2d &) const;
- };
-
-
- /**
- * This is a helper function for the write_gmv() function. There,
- * the data in the patches needs to be copied around as output is
- * one variable globally at a time, rather than all data on each
- * vertex at a time. This copying around can be done detached from
- * the main thread, and is thus moved into this separate function.
- *
- * Note that because of the similarity of the formats, this function
- * is also used by the Vtk and Tecplot output functions.
- */
- template
- static void
- write_gmv_reorder_data_vectors (const std::vector > &patches,
- Table<2,double> &data_vectors);
-
-
-
-};
+}
/**
- * This class is the interface to the DataOutBase class, as already its name
+ * This class is the interface to the functions in the DataOutBase namespace, as already its name
* might suggest. It does not offer much functionality apart from a way to
* access the implemented formats and a way to dynamically dispatch what output
* format to chose.
@@ -2278,36 +2138,9 @@ private:
* @author Wolfgang Bangerth, 1999
*/
template
-class DataOutInterface : private DataOutBase
+class DataOutInterface
{
public:
- /*
- * Import a few names that were
- * previously in this class and have then
- * moved to the base class. Since the
- * base class is inherited from
- * privately, we need to re-import these
- * symbols to make sure that references
- * to DataOutInterface::XXX
- * remain valid.
- */
- using DataOutBase::OutputFormat;
- using DataOutBase::default_format;
- using DataOutBase::dx;
- using DataOutBase::gnuplot;
- using DataOutBase::povray;
- using DataOutBase::eps;
- using DataOutBase::tecplot;
- using DataOutBase::tecplot_binary;
- using DataOutBase::vtk;
- using DataOutBase::vtu;
- using DataOutBase::deal_II_intermediate;
- using DataOutBase::parse_output_format;
- using DataOutBase::get_output_format_names;
- using DataOutBase::determine_intermediate_format_dimensions;
-
- using DataOutBase::Patch;
-
/**
* Constructor.
*/
@@ -2630,7 +2463,7 @@ public:
* the mesh and solution data were written to a single file. See
* write_xdmf_file() for an example of usage.
*/
- XDMFEntry create_xdmf_entry (const DataOutFilter &data_filter,
+ XDMFEntry create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter,
const std::string &h5_filename,
const double cur_time,
MPI_Comm comm) const;
@@ -2640,7 +2473,7 @@ public:
* the mesh and solution data were written to separate files. See
* write_xdmf_file() for an example of usage.
*/
- XDMFEntry create_xdmf_entry (const DataOutFilter &data_filter,
+ XDMFEntry create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter,
const std::string &h5_mesh_filename,
const std::string &h5_solution_filename,
const double cur_time,
@@ -2674,7 +2507,8 @@ public:
* single HDF5 file containing both the mesh and solution values.
@deprecated: use write_hdf5_parallel(DataOutFilter, ...) instead
*/
- void write_hdf5_parallel (const std::string &filename, MPI_Comm comm) const DEAL_II_DEPRECATED;
+ void write_hdf5_parallel (const std::string &filename,
+ MPI_Comm comm) const DEAL_II_DEPRECATED;
/**
* Write the data in data_filter to a single HDF5 file containing both
@@ -2689,7 +2523,8 @@ public:
* data_out.write_hdf5_parallel(data_filter, "solution.h5", MPI_COMM_WORLD);
* @endcode
*/
- void write_hdf5_parallel (const DataOutFilter &data_filter, const std::string &filename, MPI_Comm comm) const;
+ void write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter,
+ const std::string &filename, MPI_Comm comm) const;
/**
* Write the data in data_filter to HDF5 file(s). If write_mesh_file
@@ -2698,7 +2533,8 @@ public:
* is true and the filenames are the same, the resulting file will
* contain both mesh data and solution values.
*/
- void write_hdf5_parallel (const DataOutFilter &data_filter, const bool write_mesh_file, const std::string &mesh_filename, const std::string &solution_filename, MPI_Comm comm) const;
+ void write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter,
+ const bool write_mesh_file, const std::string &mesh_filename, const std::string &solution_filename, MPI_Comm comm) const;
/**
* DataOutFilter is an intermediate data format that reduces the amount of
@@ -2706,7 +2542,7 @@ public:
* can then later be used again to write data in a concrete file format;
* see, for example, DataOutBase::write_hdf5_parallel().
*/
- void write_filtered_data (DataOutFilter &filtered_data) const;
+ void write_filtered_data (DataOutBase::DataOutFilter &filtered_data) const;
/**
@@ -2724,64 +2560,64 @@ public:
* format is default_format.
*/
void write (std::ostream &out,
- const OutputFormat output_format = default_format) const;
+ const DataOutBase::OutputFormat output_format = DataOutBase::default_format) const;
/**
* Set the default format. The value set here is used anytime,
* output for format default_format is requested.
*/
- void set_default_format (const OutputFormat default_format);
+ void set_default_format (const DataOutBase::OutputFormat default_format);
/**
* Set the flags to be used for output in OpenDX format.
*/
- void set_flags (const DXFlags &dx_flags);
+ void set_flags (const DataOutBase::DXFlags &dx_flags);
/**
* Set the flags to be used for output in UCD format.
*/
- void set_flags (const UcdFlags &ucd_flags);
+ void set_flags (const DataOutBase::UcdFlags &ucd_flags);
/**
* Set the flags to be used for output in GNUPLOT format.
*/
- void set_flags (const GnuplotFlags &gnuplot_flags);
+ void set_flags (const DataOutBase::GnuplotFlags &gnuplot_flags);
/**
* Set the flags to be used for output in POVRAY format.
*/
- void set_flags (const PovrayFlags &povray_flags);
+ void set_flags (const DataOutBase::PovrayFlags &povray_flags);
/**
* Set the flags to be used for output in EPS output.
*/
- void set_flags (const EpsFlags &eps_flags);
+ void set_flags (const DataOutBase::EpsFlags &eps_flags);
/**
* Set the flags to be used for output in GMV format.
*/
- void set_flags (const GmvFlags &gmv_flags);
+ void set_flags (const DataOutBase::GmvFlags &gmv_flags);
/**
* Set the flags to be used for output in Tecplot format.
*/
- void set_flags (const TecplotFlags &tecplot_flags);
+ void set_flags (const DataOutBase::TecplotFlags &tecplot_flags);
/**
* Set the flags to be used for output in VTK format.
*/
- void set_flags (const VtkFlags &vtk_flags);
+ void set_flags (const DataOutBase::VtkFlags &vtk_flags);
/**
* Set the flags to be used for output in SVG format.
*/
- void set_flags (const SvgFlags &svg_flags);
+ void set_flags (const DataOutBase::SvgFlags &svg_flags);
/**
* Set the flags to be used for output in deal.II intermediate
* format.
*/
- void set_flags (const Deal_II_IntermediateFlags &deal_II_intermediate_flags);
+ void set_flags (const DataOutBase::Deal_II_IntermediateFlags &deal_II_intermediate_flags);
/**
* A function that returns the same string as the respective
@@ -2792,7 +2628,7 @@ public:
* calling this function.
*/
std::string
- default_suffix (const OutputFormat output_format = default_format) const;
+ default_suffix (const DataOutBase::OutputFormat output_format = DataOutBase::default_format) const;
/**
* Declare parameters for all output formats by declaring
@@ -2840,7 +2676,7 @@ protected:
* functions to know what they shall print.
*/
virtual
- const std::vector > &
+ const std::vector > &
get_patches () const = 0;
/**
@@ -2887,73 +2723,73 @@ private:
* default_format is requested. It can be changed by the
* set_format function or in a parameter file.
*/
- OutputFormat default_fmt;
+ DataOutBase::OutputFormat default_fmt;
/**
* Flags to be used upon output of OpenDX data. Can be changed by
* using the set_flags function.
*/
- DXFlags dx_flags;
+ DataOutBase::DXFlags dx_flags;
/**
* Flags to be used upon output of UCD data. Can be changed by using
* the set_flags function.
*/
- UcdFlags ucd_flags;
+ DataOutBase::UcdFlags ucd_flags;
/**
* Flags to be used upon output of GNUPLOT data. Can be changed by
* using the set_flags function.
*/
- GnuplotFlags gnuplot_flags;
+ DataOutBase::GnuplotFlags gnuplot_flags;
/**
* Flags to be used upon output of POVRAY data. Can be changed by
* using the set_flags function.
*/
- PovrayFlags povray_flags;
+ DataOutBase::PovrayFlags povray_flags;
/**
* Flags to be used upon output of EPS data in one space
* dimension. Can be changed by using the set_flags
* function.
*/
- EpsFlags eps_flags;
+ DataOutBase::EpsFlags eps_flags;
/**
* Flags to be used upon output of gmv data in one space
* dimension. Can be changed by using the set_flags
* function.
*/
- GmvFlags gmv_flags;
+ DataOutBase::GmvFlags gmv_flags;
/**
* Flags to be used upon output of Tecplot data in one space
* dimension. Can be changed by using the set_flags
* function.
*/
- TecplotFlags tecplot_flags;
+ DataOutBase::TecplotFlags tecplot_flags;
/**
* Flags to be used upon output of vtk data in one space
* dimension. Can be changed by using the set_flags
* function.
*/
- VtkFlags vtk_flags;
+ DataOutBase::VtkFlags vtk_flags;
/**
* Flags to be used upon output of svg data in one space
* dimension. Can be changed by using the set_flags
* function.
*/
- SvgFlags svg_flags;
+ DataOutBase::SvgFlags svg_flags;
/**
* Flags to be used upon output of deal.II intermediate data in one
* space dimension. Can be changed by using the set_flags
* function.
*/
- Deal_II_IntermediateFlags deal_II_intermediate_flags;
+ DataOutBase::Deal_II_IntermediateFlags deal_II_intermediate_flags;
};
@@ -3077,7 +2913,7 @@ protected:
* It returns the patches as read the last time a stream was given
* to the read() function.
*/
- virtual const std::vector > &
+ virtual const std::vector > &
get_patches () const;
/**
@@ -3117,7 +2953,7 @@ private:
* Arrays holding the set of patches as well as the names of output
* variables, all of which we read from an input stream.
*/
- std::vector > patches;
+ std::vector > patches;
std::vector dataset_names;
/**
diff --git a/deal.II/source/base/data_out_base.cc b/deal.II/source/base/data_out_base.cc
index 23861d2c09..ef6471dbd3 100644
--- a/deal.II/source/base/data_out_base.cc
+++ b/deal.II/source/base/data_out_base.cc
@@ -273,6 +273,146 @@ namespace
#endif
}
+
+// some declarations of functions
+namespace DataOutBase
+{
+ /**
+ * Write the coordinates of nodes in the desired format.
+ */
+ template
+ static void write_nodes (const std::vector > &patches,
+ STREAM &out);
+
+ /**
+ * Write the node numbers of a cell in the desired format.
+ */
+ template
+ static void write_cells (const std::vector > &patches,
+ STREAM &out);
+
+ /**
+ * Write data in the desired format.
+ */
+ template
+ static void write_data (const std::vector > &patches,
+ const unsigned int n_data_sets,
+ const bool double_precision,
+ STREAM &out);
+
+
+ /**
+ * This function projects a three-dimensional point (Point<3> point)
+ * onto a two-dimensional image plane, specified by the position of
+ * the camera viewing system (Point<3> camera_position), camera
+ * direction (Point<3> camera_position), camera horizontal (Point<3>
+ * camera_horizontal, necessary for the correct alignment of the
+ * later images), and the focus of the camera (float camera_focus).
+ *
+ * For SVG output.
+ */
+ static Point<2> svg_project_point(Point<3> point,
+ Point<3> camera_position,
+ Point<3> camera_direction,
+ Point<3> camera_horizontal,
+ float camera_focus);
+ /**
+ * Function to compute the gradient parameters for a triangle with
+ * given values for the vertices.
+ *
+ * Used for svg output.
+ */
+ static Point<6> svg_get_gradient_parameters(Point<3> points[]);
+
+ /**
+ * Class holding the data of one cell of a patch in two space
+ * dimensions for output. It is the projection of a cell in
+ * three-dimensional space (two coordinates, one height value) to
+ * the direction of sight.
+ */
+ class SvgCell
+ {
+ public:
+
+ // Center of the cell (three-dimensional)
+ Point<3> center;
+
+ /**
+ * Vector of vertices of this cell (three-dimensional)
+ */
+ Point<3> vertices[4];
+
+ /**
+ * Depth into the picture, which is defined as the distance from
+ * an observer at an the origin in direction of the line of sight.
+ */
+ float depth;
+
+ /**
+ * Vector of vertices of this cell (projected, two-dimensional).
+ */
+ Point<2> projected_vertices[4];
+
+ // Center of the cell (projected, two-dimensional)
+ Point<2> projected_center;
+
+ /**
+ * Comparison operator for sorting.
+ */
+ bool operator < (const SvgCell &) const;
+ };
+
+
+ /**
+ * Class holding the data of one cell of a patch in two space
+ * dimensions for output. It is the projection of a cell in
+ * three-dimensional space (two coordinates, one height value) to
+ * the direction of sight.
+ */
+ class EpsCell2d
+ {
+ public:
+
+ /**
+ * Vector of vertices of this cell.
+ */
+ Point<2> vertices[4];
+
+ /**
+ * Data value from which the actual colors will be computed by the
+ * colorization function stated in the EpsFlags class.
+ */
+ float color_value;
+
+ /**
+ * Depth into the picture, which is defined as the distance from
+ * an observer at an the origin in direction of the line of sight.
+ */
+ float depth;
+
+ /**
+ * Comparison operator for sorting.
+ */
+ bool operator < (const EpsCell2d &) const;
+ };
+
+
+ /**
+ * This is a helper function for the write_gmv() function. There,
+ * the data in the patches needs to be copied around as output is
+ * one variable globally at a time, rather than all data on each
+ * vertex at a time. This copying around can be done detached from
+ * the main thread, and is thus moved into this separate function.
+ *
+ * Note that because of the similarity of the formats, this function
+ * is also used by the Vtk and Tecplot output functions.
+ */
+ template
+ static void
+ write_gmv_reorder_data_vectors (const std::vector > &patches,
+ Table<2,double> &data_vectors);
+}
+
//----------------------------------------------------------------------//
// DataOutFilter class member functions
//----------------------------------------------------------------------//
@@ -5629,16 +5769,6 @@ void DataOutBase::write_vtu_main (const std::vector > &patch
AssertThrow (out, ExcIO());
}
-template <>
-void DataOutBase::write_svg<1,1> (const std::vector > &patches,
- const std::vector &data_names,
- const std::vector > &vector_data_ranges,
- const SvgFlags &flags,
- std::ostream &out)
-{
- AssertThrow (false, ExcNotImplemented());
-}
-
template
void DataOutBase::write_svg (const std::vector > &patches,
@@ -5971,7 +6101,7 @@ void DataOutBase::write_svg (const std::vector > &patches,
}
-// write the svg file
+ // write the svg file
if (width==0)
width = static_cast(.5 + height * (x_dimension_perspective / y_dimension_perspective));
unsigned int additional_width = 0;
@@ -6780,16 +6910,25 @@ create_xdmf_entry (const std::string &h5_filename, const double cur_time, MPI_Co
return create_xdmf_entry(data_filter, h5_filename, cur_time, comm);
}
+
+
template
XDMFEntry DataOutInterface::
-create_xdmf_entry (const DataOutFilter &data_filter, const std::string &h5_filename, const double cur_time, MPI_Comm comm) const
+create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter,
+ const std::string &h5_filename, const double cur_time, MPI_Comm comm) const
{
return create_xdmf_entry(data_filter, h5_filename, h5_filename, cur_time, comm);
}
+
+
template
XDMFEntry DataOutInterface::
-create_xdmf_entry (const DataOutFilter &data_filter, const std::string &h5_mesh_filename, const std::string &h5_solution_filename, const double cur_time, MPI_Comm comm) const
+create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter,
+ const std::string &h5_mesh_filename,
+ const std::string &h5_solution_filename,
+ const double cur_time,
+ MPI_Comm comm) const
{
unsigned int local_node_cell_count[2], global_node_cell_count[2];
int myrank;
@@ -6936,7 +7075,7 @@ std::string XDMFEntry::get_xdmf_content(const unsigned int indent_level) const
*/
template
void DataOutInterface::
-write_filtered_data (DataOutFilter &filtered_data) const
+write_filtered_data (DataOutBase::DataOutFilter &filtered_data) const
{
DataOutBase::write_filtered_data(get_patches(), get_dataset_names(),
get_vector_data_ranges(),
@@ -6947,7 +7086,7 @@ template
void DataOutBase::write_filtered_data (const std::vector > &patches,
const std::vector &data_names,
const std::vector > &vector_data_ranges,
- DataOutFilter &filtered_data)
+ DataOutBase::DataOutFilter &filtered_data)
{
const unsigned int n_data_sets = data_names.size();
unsigned int n_node, n_cell;
@@ -7050,21 +7189,23 @@ write_hdf5_parallel (const std::string &filename, MPI_Comm comm) const
template
void DataOutInterface::
-write_hdf5_parallel (const DataOutFilter &data_filter, const std::string &filename, MPI_Comm comm) const
+write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter,
+ const std::string &filename, MPI_Comm comm) const
{
DataOutBase::write_hdf5_parallel(get_patches(), data_filter, filename, comm);
}
template
void DataOutInterface::
-write_hdf5_parallel (const DataOutFilter &data_filter, const bool write_mesh_file, const std::string &mesh_filename, const std::string &solution_filename, MPI_Comm comm) const
+write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter,
+ const bool write_mesh_file, const std::string &mesh_filename, const std::string &solution_filename, MPI_Comm comm) const
{
DataOutBase::write_hdf5_parallel(get_patches(), data_filter, write_mesh_file, mesh_filename, solution_filename, comm);
}
template
void DataOutBase::write_hdf5_parallel (const std::vector > &patches,
- const DataOutFilter &data_filter,
+ const DataOutBase::DataOutFilter &data_filter,
const std::string &filename,
MPI_Comm comm)
{
@@ -7073,7 +7214,7 @@ void DataOutBase::write_hdf5_parallel (const std::vector > &
template
void DataOutBase::write_hdf5_parallel (const std::vector > &patches,
- const DataOutFilter &data_filter,
+ const DataOutBase::DataOutFilter &data_filter,
const bool write_mesh_file,
const std::string &mesh_filename,
const std::string &solution_filename,
@@ -7356,62 +7497,62 @@ void DataOutBase::write_hdf5_parallel (const std::vector > &
template
void
DataOutInterface::write (std::ostream &out,
- const OutputFormat output_format_) const
+ const DataOutBase::OutputFormat output_format_) const
{
- OutputFormat output_format = output_format_;
- if (output_format == default_format)
+ DataOutBase::OutputFormat output_format = output_format_;
+ if (output_format == DataOutBase::default_format)
output_format = default_fmt;
switch (output_format)
{
- case none:
+ case DataOutBase::none:
break;
- case dx:
+ case DataOutBase::dx:
write_dx (out);
break;
- case ucd:
+ case DataOutBase::ucd:
write_ucd (out);
break;
- case gnuplot:
+ case DataOutBase::gnuplot:
write_gnuplot (out);
break;
- case povray:
+ case DataOutBase::povray:
write_povray (out);
break;
- case eps:
+ case DataOutBase::eps:
write_eps (out);
break;
- case gmv:
+ case DataOutBase::gmv:
write_gmv (out);
break;
- case tecplot:
+ case DataOutBase::tecplot:
write_tecplot (out);
break;
- case tecplot_binary:
+ case DataOutBase::tecplot_binary:
write_tecplot_binary (out);
break;
- case vtk:
+ case DataOutBase::vtk:
write_vtk (out);
break;
- case vtu:
+ case DataOutBase::vtu:
write_vtu (out);
break;
- case svg:
+ case DataOutBase::svg:
write_svg (out);
break;
- case deal_II_intermediate:
+ case DataOutBase::deal_II_intermediate:
write_deal_II_intermediate (out);
break;
@@ -7424,9 +7565,9 @@ DataOutInterface::write (std::ostream &out,
template
void
-DataOutInterface::set_default_format(const OutputFormat fmt)
+DataOutInterface::set_default_format(const DataOutBase::OutputFormat fmt)
{
- Assert (fmt != default_format, ExcNotImplemented());
+ Assert (fmt != DataOutBase::default_format, ExcNotImplemented());
default_fmt = fmt;
}
@@ -7434,7 +7575,7 @@ DataOutInterface::set_default_format(const OutputFormat fmt)
template
void
-DataOutInterface::set_flags (const DXFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::DXFlags &flags)
{
dx_flags = flags;
}
@@ -7443,7 +7584,7 @@ DataOutInterface::set_flags (const DXFlags &flags)
template
void
-DataOutInterface::set_flags (const UcdFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::UcdFlags &flags)
{
ucd_flags = flags;
}
@@ -7452,7 +7593,7 @@ DataOutInterface::set_flags (const UcdFlags &flags)
template
void
-DataOutInterface::set_flags (const GnuplotFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::GnuplotFlags &flags)
{
gnuplot_flags = flags;
}
@@ -7461,7 +7602,7 @@ DataOutInterface::set_flags (const GnuplotFlags &flags)
template
void
-DataOutInterface::set_flags (const PovrayFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::PovrayFlags &flags)
{
povray_flags = flags;
}
@@ -7470,7 +7611,7 @@ DataOutInterface::set_flags (const PovrayFlags &flags)
template
void
-DataOutInterface::set_flags (const EpsFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::EpsFlags &flags)
{
eps_flags = flags;
}
@@ -7479,7 +7620,7 @@ DataOutInterface::set_flags (const EpsFlags &flags)
template
void
-DataOutInterface::set_flags (const GmvFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::GmvFlags &flags)
{
gmv_flags = flags;
}
@@ -7488,7 +7629,7 @@ DataOutInterface::set_flags (const GmvFlags &flags)
template
void
-DataOutInterface::set_flags (const TecplotFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::TecplotFlags &flags)
{
tecplot_flags = flags;
}
@@ -7497,7 +7638,7 @@ DataOutInterface::set_flags (const TecplotFlags &flags)
template
void
-DataOutInterface::set_flags (const VtkFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::VtkFlags &flags)
{
vtk_flags = flags;
}
@@ -7506,7 +7647,7 @@ DataOutInterface::set_flags (const VtkFlags &flags)
template
void
-DataOutInterface::set_flags (const SvgFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::SvgFlags &flags)
{
svg_flags = flags;
}
@@ -7515,7 +7656,7 @@ DataOutInterface::set_flags (const SvgFlags &flags)
template
void
-DataOutInterface::set_flags (const Deal_II_IntermediateFlags &flags)
+DataOutInterface::set_flags (const DataOutBase::Deal_II_IntermediateFlags &flags)
{
deal_II_intermediate_flags = flags;
}
@@ -7525,9 +7666,9 @@ DataOutInterface::set_flags (const Deal_II_IntermediateFlags &flag
template
std::string
DataOutInterface::
-default_suffix (const OutputFormat output_format) const
+default_suffix (const DataOutBase::OutputFormat output_format) const
{
- if (output_format == default_format)
+ if (output_format == DataOutBase::default_format)
return DataOutBase::default_suffix (default_fmt);
else
return DataOutBase::default_suffix (output_format);
@@ -7540,46 +7681,46 @@ void
DataOutInterface::declare_parameters (ParameterHandler &prm)
{
prm.declare_entry ("Output format", "gnuplot",
- Patterns::Selection (get_output_format_names ()),
+ Patterns::Selection (DataOutBase::get_output_format_names ()),
"A name for the output format to be used");
prm.declare_entry("Subdivisions", "1", Patterns::Integer(),
"Number of subdivisions of each mesh cell");
prm.enter_subsection ("DX output parameters");
- DXFlags::declare_parameters (prm);
+ DataOutBase::DXFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("UCD output parameters");
- UcdFlags::declare_parameters (prm);
+ DataOutBase::UcdFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("Gnuplot output parameters");
- GnuplotFlags::declare_parameters (prm);
+ DataOutBase::GnuplotFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("Povray output parameters");
- PovrayFlags::declare_parameters (prm);
+ DataOutBase::PovrayFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("Eps output parameters");
- EpsFlags::declare_parameters (prm);
+ DataOutBase::EpsFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("Gmv output parameters");
- GmvFlags::declare_parameters (prm);
+ DataOutBase::GmvFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("Tecplot output parameters");
- TecplotFlags::declare_parameters (prm);
+ DataOutBase::TecplotFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("Vtk output parameters");
- VtkFlags::declare_parameters (prm);
+ DataOutBase::VtkFlags::declare_parameters (prm);
prm.leave_subsection ();
prm.enter_subsection ("deal.II intermediate output parameters");
- Deal_II_IntermediateFlags::declare_parameters (prm);
+ DataOutBase::Deal_II_IntermediateFlags::declare_parameters (prm);
prm.leave_subsection ();
}
@@ -7590,7 +7731,7 @@ void
DataOutInterface::parse_parameters (ParameterHandler &prm)
{
const std::string &output_name = prm.get ("Output format");
- default_fmt = parse_output_format (output_name);
+ default_fmt = DataOutBase::parse_output_format (output_name);
default_subdivisions = prm.get_integer ("Subdivisions");
prm.enter_subsection ("DX output parameters");
@@ -7695,7 +7836,7 @@ DataOutReader::read (std::istream &in)
{
std::pair
dimension_info
- = this->determine_intermediate_format_dimensions (in);
+ = DataOutBase::determine_intermediate_format_dimensions (in);
AssertThrow ((dimension_info.first == dim) &&
(dimension_info.second == spacedim),
ExcIncompatibleDimensions (dimension_info.first, dim,