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 | Private Attributes | List of all members
DataOutReader< dim, spacedim > Class Template Reference

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

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

Public Member Functions

void read (std::istream &in)
 
void merge (const DataOutReader< dim, spacedim > &other)
 
- Public Member Functions inherited from DataOutInterface< dim, dim >
 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)
 
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 ::ExceptionBaseExcIncompatibleDatasetNames ()
 
static ::ExceptionBaseExcIncompatiblePatchLists ()
 
static ::ExceptionBaseExcIncompatibleDimensions (int arg1, int arg2, int arg3, int arg4)
 
- Static Public Member Functions inherited from DataOutInterface< dim, dim >
static void declare_parameters (ParameterHandler &prm)
 

Protected Member Functions

virtual const std::vector<::DataOutBase::Patch< dim, spacedim > > & get_patches () const override
 
virtual std::vector< std::string > get_dataset_names () const override
 
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > get_nonscalar_data_ranges () const override
 
- Protected Member Functions inherited from DataOutInterface< dim, dim >
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
 

Private Attributes

std::vector<::DataOutBase::Patch< dim, spacedim > > patches
 
std::vector< std::string > dataset_names
 
std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > nonscalar_data_ranges
 

Additional Inherited Members

- Protected Attributes inherited from DataOutInterface< dim, dim >
unsigned int default_subdivisions
 

Detailed Description

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

A class that is used to read data written in deal.II intermediate format back in, so that it can be written out in any of the other supported graphics formats. This class has two main purposes:

The first use of this class is so that application programs can defer the decision of which graphics format to use until after the program has been run. The data is written in intermediate format into a file, and later on it can then be converted into any graphics format you wish. This may be useful, for example, if you want to convert it to gnuplot format to get a quick glimpse and later on want to convert it to OpenDX format as well to get a high quality version of the data. The present class allows to read this intermediate format back into the program, and allows it to be written in any other supported format using the relevant functions of the base class.

The second use is mostly useful in parallel programs: rather than having one central process generate the graphical output for the entire program, one can let each process generate the graphical data for the cells it owns, and write it into a separate file in intermediate format. Later on, all these intermediate files can then be read back in and merged together, a process that is fast compared to generating the data in the first place. The use of the intermediate format is mostly because it allows separate files to be merged, while this is almost impossible once the data has been written out in any of the supported established graphics formats.

This second use scenario is explained in some detail in the step-18 example program.

In order to read data back into this object, you have to know the template parameters for the space dimension which were used when writing the data. If this knowledge is available at compile time, then this is no problem. However, if it is not (such as in a simple format converter), then it needs to be figured out at run time, even though the compiler already needs it at compile time. A way around using the DataOutBase::determine_intermediate_format_dimensions() function.

Note that the intermediate format is what its name suggests: a direct representation of internal data. It isn't standardized and will change whenever we change our internal representation. You can only expect to process files written in this format using the same version of deal.II that was used for writing.

Author
Wolfgang Bangerth, 2005

Definition at line 3188 of file data_out_base.h.


The documentation for this class was generated from the following files: