Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
DataOutInterface< dim, spacedim > Class Template Referenceabstract

#include <deal.II/base/data_out_base.h>

Inheritance diagram for DataOutInterface< dim, spacedim >:
[legend]

Public Member Functions

 DataOutInterface ()
 
virtual ~DataOutInterface ()=default
 
void write_dx (std::ostream &out) const
 
void write_eps (std::ostream &out) const
 
void write_gmv (std::ostream &out) const
 
void write_gnuplot (std::ostream &out) const
 
void write_povray (std::ostream &out) const
 
void write_tecplot (std::ostream &out) const
 
void write_tecplot_binary (std::ostream &out) const
 
void write_ucd (std::ostream &out) const
 
void write_vtk (std::ostream &out) const
 
void write_vtu (std::ostream &out) const
 
void write_vtu_in_parallel (const std::string &filename, MPI_Comm comm) const
 
void write_pvtu_record (std::ostream &out, const std::vector< std::string > &piece_names) const
 
std::string write_vtu_with_pvtu_record (const std::string &directory, const std::string &filename_without_extension, const unsigned int counter, const MPI_Comm &mpi_communicator, const unsigned int n_digits_for_counter=numbers::invalid_unsigned_int, const unsigned int n_groups=0) const
 
void write_svg (std::ostream &out) const
 
void write_deal_II_intermediate (std::ostream &out) const
 
XDMFEntry create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter, const std::string &h5_filename, const double cur_time, MPI_Comm comm) const
 
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, MPI_Comm comm) const
 
void write_xdmf_file (const std::vector< XDMFEntry > &entries, 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
 
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
 
void write_filtered_data (DataOutBase::DataOutFilter &filtered_data) const
 
void write (std::ostream &out, const DataOutBase::OutputFormat output_format=DataOutBase::default_format) const
 
void set_default_format (const DataOutBase::OutputFormat default_format)
 
template<typename FlagType >
void set_flags (const FlagType &flags)
 
std::string default_suffix (const DataOutBase::OutputFormat output_format=DataOutBase::default_format) const
 
void parse_parameters (ParameterHandler &prm)
 
std::size_t memory_consumption () const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Protected Member Functions

virtual const std::vector< DataOutBase::Patch< dim, spacedim > > & get_patches () const =0
 
virtual std::vector< std::string > get_dataset_names () const =0
 
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > get_nonscalar_data_ranges () const
 
void validate_dataset_names () const
 

Protected Attributes

unsigned int default_subdivisions
 

Private Attributes

DataOutBase::OutputFormat default_fmt
 
DataOutBase::DXFlags dx_flags
 
DataOutBase::UcdFlags ucd_flags
 
DataOutBase::GnuplotFlags gnuplot_flags
 
DataOutBase::PovrayFlags povray_flags
 
DataOutBase::EpsFlags eps_flags
 
DataOutBase::GmvFlags gmv_flags
 
DataOutBase::TecplotFlags tecplot_flags
 
DataOutBase::VtkFlags vtk_flags
 
DataOutBase::SvgFlags svg_flags
 
DataOutBase::Deal_II_IntermediateFlags deal_II_intermediate_flags
 

Detailed Description

template<int dim, int spacedim = dim>
class DataOutInterface< dim, spacedim >

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.

This class is thought as a base class to classes actually generating data for output. It has two abstract virtual functions, get_patches() and get_dataset_names() produce the data which is actually needed. These are the only functions that need to be overloaded by a derived class. In additional to that, it has a function for each output format supported by the underlying base class which gets the output data using these two virtual functions and passes them to the raw output functions.

The purpose of this class is mainly two-fold: to support storing flags by which the output in the different output formats are controlled, and means to work with output in a way where output format, flags and other things are determined at run time. In addition to that it offers the abstract interface to derived classes briefly discussed above.

Output flags

The way we treat flags in this class is very similar to that used in the GridOut class. For detailed information on the why's and how's, as well as an example of programming, we refer to the documentation of that class.

Basically, this class stores a set of flags for each output format supported by the underlying DataOutBase class. These are used whenever one of the write_* functions is used. By default, the values of these flags are set to reasonable start-ups, but in case you want to change them, you can create a structure holding the flags for one of the output formats and set it using the set_flags functions of this class to determine all future output the object might produce by that output format.

For information on what parameters are supported by different output functions, please see the documentation of the DataOutBase class and its member classes.

Run time selection of output parameters

In the output flags classes, described above, many flags are defined for output in the different formats. In order to make them available to the input file handler class ParameterHandler, each of these has a function declaring these flags to the parameter handler and to read them back from an actual input file. In order to avoid that in user programs these functions have to be called for each available output format and the respective flag class, the present DataOutInterface class offers a function declare_parameters which calls the respective function of all known output format flags classes. The flags of each such format are packed together in a subsection in the input file. Likewise, there is a function parse_parameters which reads these parameters and stores them in the flags associated with this object (see above).

Using these functions, you do not have to track which formats are presently implemented.

Usage is as follows:

// within function declaring parameters:
prm.enter_subsection("Output format options");
prm.leave_subsection();
...
// within function doing the output:
DataOut<dim> out;
prm.enter_subsection("Output format options");
out.parse_parameters(prm);
prm.leave_subsection();

Note that in the present example, the class DataOut was used. However, any other class derived from DataOutInterface would work alike.

Run time selection of formats

This class, much like the GridOut class, has a set of functions providing a list of supported output formats, an enum denoting all these and a function to parse a string and return the respective enum value if it is a valid output format's name (actually, these functions are inherited from the base class). Finally, there is a function write, which takes a value of this enum and dispatches to one of the actual write_* functions depending on the output format selected by this value.

The functions offering the different output format names are, respectively, default_suffix, parse_output_format, and get_output_format_names. They make the selection of output formats in parameter files much easier, and especially independent of the formats presently implemented. User programs need therefore not be changed whenever a new format is implemented.

Additionally, objects of this class have a default format, which can be set by the parameter "Output format" of the parameter file. Within a program, this can be changed by the member function set_default_format. Using this default format, it is possible to leave the format selection completely to the parameter file. A suitable suffix for the output file name can be obtained by default_suffix without arguments.

Author
Wolfgang Bangerth, 1999, Denis Davydov, 2018

Definition at line 2598 of file data_out_base.h.


The documentation for this class was generated from the following files:
DataOutInterface::declare_parameters
static void declare_parameters(ParameterHandler &prm)
Definition: data_out_base.cc:7885