Reference documentation for deal.II version 9.3.3
\(\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 Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
DataOut< dim, DoFHandlerType > Class Template Reference

#include <deal.II/numerics/data_out.h>

Inheritance diagram for DataOut< dim, DoFHandlerType >:
[legend]

Public Types

enum  CurvedCellRegion { no_curved_cells , curved_boundary , curved_inner_cells }
 
using cell_iterator = typename DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension >::cell_iterator
 
using FirstCellFunctionType = typename std::function< cell_iterator(const Triangulation< dim, spacedim > &)>
 
using NextCellFunctionType = typename std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)>
 
enum  DataVectorType
 

Public Member Functions

 DataOut ()
 
virtual void build_patches (const unsigned int n_subdivisions=0)
 
virtual void build_patches (const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const unsigned int n_subdivisions=0, const CurvedCellRegion curved_region=curved_boundary)
 
virtual void build_patches (const hp::MappingCollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const unsigned int n_subdivisions=0, const CurvedCellRegion curved_region=curved_boundary)
 
void set_cell_selection (const std::function< cell_iterator(const Triangulation< dim, spacedim > &)> &first_cell, const std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> &next_cell)
 
void set_cell_selection (const FilteredIterator< cell_iterator > &filtered_iterator)
 
const std::pair< FirstCellFunctionType, NextCellFunctionTypeget_cell_selection () const
 
virtual cell_iterator first_cell ()
 
virtual cell_iterator next_cell (const cell_iterator &cell)
 
void attach_dof_handler (const DoFHandler< dim > &)
 
void attach_triangulation (const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &)
 
void add_data_vector (const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_data_vector (const VectorType &data, const std::string &name, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_data_vector (const DoFHandler< dim > &dof_handler, const VectorType &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_data_vector (const DoFHandler< dim > &dof_handler, const VectorType &data, const std::string &name, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_data_vector (const VectorType &data, const DataPostprocessor< DoFHandlerType::space_dimension > &data_postprocessor)
 
void add_data_vector (const DoFHandler< dim > &dof_handler, const VectorType &data, const DataPostprocessor< DoFHandlerType::space_dimension > &data_postprocessor)
 
void add_mg_data_vector (const DoFHandler< dim > &dof_handler, const MGLevelObject< VectorType > &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_mg_data_vector (const DoFHandler< dim > &dof_handler, const MGLevelObject< VectorType > &data, const std::string &name)
 
void clear_data_vectors ()
 
void clear_input_data_references ()
 
void merge_patches (const DataOut_DoFData< DoFHandlerType2, patch_dim, patch_space_dim > &source, const Point< patch_space_dim > &shift=Point< patch_space_dim >())
 
virtual void clear ()
 
std::size_t memory_consumption () const
 
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_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, const 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, const 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, const MPI_Comm &comm) const
 
void write_xdmf_file (const std::vector< XDMFEntry > &entries, const std::string &filename, const MPI_Comm &comm) const
 
void write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter, const std::string &filename, const 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, const 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)
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Static Public Attributes

static constexpr unsigned int spacedim = DoFHandlerType::space_dimension
 

Protected Types

using Patch = ::DataOutBase::Patch< patch_dim, patch_space_dim >
 

Protected Member Functions

virtual const std::vector< Patch > & get_patches () const override
 
virtual std::vector< std::string > get_dataset_names () const override
 
std::vector< std::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > get_fes () const
 
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > get_nonscalar_data_ranges () const override
 
void validate_dataset_names () const
 

Protected Attributes

SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > triangulation
 
SmartPointer< const DoFHandler< dim > > dofs
 
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > dof_data
 
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > cell_data
 
std::vector< Patchpatches
 
unsigned int default_subdivisions
 

Private Member Functions

virtual cell_iterator first_locally_owned_cell ()
 
virtual cell_iterator next_locally_owned_cell (const cell_iterator &cell)
 
void build_one_patch (const std::pair< cell_iterator, unsigned int > *cell_and_index, internal::DataOutImplementation::ParallelData< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &scratch_data, const unsigned int n_subdivisions, const CurvedCellRegion curved_cell_region)
 
void add_data_vector_internal (const DoFHandler< dim > *dof_handler, const VectorType &data, const std::vector< std::string > &names, const DataVectorType type, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation, const bool deduce_output_names)
 

Private Attributes

std::function< cell_iterator(const Triangulation< dim, spacedim > &)> first_cell_function
 
std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> next_cell_function
 
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, typename DoFHandlerType = DoFHandler<dim>>
class DataOut< dim, DoFHandlerType >

This class is the main class to provide output of data described by finite element fields defined on a collection of cells.

This class is an actual implementation of the functionality proposed by the DataOut_DoFData class. It offers a function build_patches() that generates the data to be written in some graphics format. Most of the interface and an example of its use is described in the documentation of this base class.

The only thing this class offers is the function build_patches() which loops over all cells of the triangulation stored by the attach_dof_handler() function of the base class (with the exception of cells of parallel::distributed::Triangulation objects that are not owned by the current processor) and converts the data on these to actual patches which are the objects that are later output by the functions of the base classes. You can give a parameter to the function which determines how many subdivisions in each coordinate direction are to be performed, i.e. of how many subcells each patch shall consist. The default is one, but you may want to choose a higher number for higher order elements, for example two for quadratic elements, three for cubic elements, and so on. (See step-11 for an example.) The purpose of this parameter is because most graphics programs do not allow to specify higher order polynomial functions in the file formats: only data at vertices can be plotted and is then shown as a bilinear interpolation within the interior of cells. This may be insufficient if you have higher order finite elements, and the only way to achieve better output is to subdivide each cell of the mesh into several cells for graphical output. Of course, what you get to see is still a bilinear interpolation on each cell of the output (where these cells are not subdivisions of the cells of the triangulation in use) due to the same limitations in output formats, but at least a bilinear interpolation of a higher order polynomial on a finer mesh.

Note that after having called build_patches() once, you can call one or more of the write() functions of DataOutInterface. You can therefore output the same data in more than one format without having to rebuild the patches.

User interface information

The base classes of this class, DataOutBase, DataOutInterface and DataOut_DoFData offer several interfaces of their own. Refer to the DataOutBase class's documentation for a discussion of the different output formats presently supported, DataOutInterface for ways of selecting which format to use upon output at run-time and without the need to adapt your program when new formats become available, as well as for flags to determine aspects of output. The DataOut_DoFData() class's documentation has an example of using nodal data to generate output.

Extensions

By default, this class produces patches for all active cells. Sometimes, this is not what you want, maybe because there are simply too many (and too small to be seen individually) or because you only want to see a certain region of the domain (for example only in the fluid part of the domain in step-46), or for some other reason.

For this, internally build_patches() does not generate the sequence of cells to be converted into patches itself, but relies on the two function that we'll call first_cell() and next_cell(). By default, they return the first active cell, and the next active cell, respectively. But this can be changed using the set_cell_selection() function that allows you to replace this behavior. What set_cell_selection() wants to know is how you want to pick out the first cell on which output should be generated, and how given one cell on which output is generated you want to pick the next cell.

This may, for example, include only cells that are in parts of a domain (e.g., if you don't care about the solution elsewhere, think for example a buffer region in which you attenuate outgoing waves in the Perfectly Matched Layer method) or if you don't want output to be generated at all levels of an adaptively refined mesh because this creates too much data (in this case, the set of cells returned by your implementations of the first_cell and next_cell arguments to set_cell_selection() will include non-active cells, and DataOut::build_patches() will simply take interpolated values of the solution instead of the exact values on these cells children for output).

Definition at line 147 of file data_out.h.

Member Typedef Documentation

◆ cell_iterator

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
using DataOut< dim, DoFHandlerType >::cell_iterator = typename DataOut_DoFData<DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension>::cell_iterator

Typedef to the iterator type of the dof handler class under consideration.

Definition at line 162 of file data_out.h.

◆ FirstCellFunctionType

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
using DataOut< dim, DoFHandlerType >::FirstCellFunctionType = typename std::function<cell_iterator(const Triangulation<dim, spacedim> &)>

The type of the function object returning the first cell as used in set_cell_selection().

Definition at line 171 of file data_out.h.

◆ NextCellFunctionType

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
using DataOut< dim, DoFHandlerType >::NextCellFunctionType = typename std::function<cell_iterator(const Triangulation<dim, spacedim> &, const cell_iterator &)>

The type of the function object returning the next cell as used in set_cell_selection().

Definition at line 178 of file data_out.h.

◆ Patch

using DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::Patch = ::DataOutBase::Patch<patch_dim, patch_space_dim>
protectedinherited

Abbreviate the somewhat lengthy name for the Patch class.

Definition at line 965 of file data_out_dof_data.h.

Member Enumeration Documentation

◆ CurvedCellRegion

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
enum DataOut::CurvedCellRegion

Enumeration describing the part of the domain in which cells should be written with curved boundaries. In reality, no file format we are aware of really supports curved boundaries, but this can be emulated by plotting edges as a sequence of straight lines (and faces in 3d as a collection of bilinear patches) if DataOut::build_patches() is called with a number of subdivisions greater than 1.

The elements of this enumeration then describe for which cells DataOut::build_patches() will query the manifold or boundary description for curved geometries.

Enumerator
no_curved_cells 

The geometry or boundary description will never be queried for curved geometries. This means that even if you have more than one subdivision per cell (see DataOut::build_patches() for what exactly this means) and even if the geometry really is curved, each cell will still be subdivided as if it was just a bi- or trilinear cell.

curved_boundary 

The geometry or boundary description will be queried for curved geometries for cells located at the boundary, i.e., for cells that have at least one face at the boundary. This is sufficient if you have not attached a manifold description to the interiors of cells but only to faces at the boundary.

curved_inner_cells 

The geometry description will be queried for all cells and all faces, whether they are at the boundary or not. This option is appropriate if you have attached a manifold object to cells (not only to boundary faces).

Definition at line 195 of file data_out.h.

◆ DataVectorType

Type describing what the vector given to add_data_vector() is: a vector that has one entry per degree of freedom in a DoFHandler object (such as solution vectors), or one entry per cell in the triangulation underlying the DoFHandler object (such as error per cell data). The value type_automatic tells add_data_vector() to find out itself (see the documentation of add_data_vector() for the method used).

Definition at line 617 of file data_out_dof_data.h.

Constructor & Destructor Documentation

◆ DataOut()

template<int dim, typename DoFHandlerType >
DataOut< dim, DoFHandlerType >::DataOut

Constructor.

Definition at line 69 of file data_out.cc.

Member Function Documentation

◆ build_patches() [1/3]

template<int dim, typename DoFHandlerType >
void DataOut< dim, DoFHandlerType >::build_patches ( const unsigned int  n_subdivisions = 0)
virtual

This is the central function of this class since it builds the list of patches to be written by the low-level functions of the base class. A patch is, in essence, some intermediate representation of the data on each cell of a triangulation and DoFHandler object that can then be used to write files in some format that is readable by visualization programs.

You can find an overview of the use of this function in the general documentation of this class. An example is also provided in the documentation of this class's base class DataOut_DoFData.

Parameters
n_subdivisionsA parameter that determines how many "patches" this function will build out of every cell. If you do not specify this value in calling, or provide the default value zero, then this is interpreted as DataOutInterface::default_subdivisions which most of the time will be equal to one (unless you have set it to something else). The purpose of this parameter is to subdivide each cell of the mesh into \(2\times 2, 3\times 3, \ldots\) "patches" in 2d, and \(2\times 2\times 2, 3\times 3\times 3, \ldots\) (if passed the value 2, 3, etc) where each patch represents the data from a regular subdivision of the cell into equal parts. Most of the times, this is not necessary and outputting one patch per cell is exactly what you want to plot the solution. That said, the data we write into files for visualization can only represent (bi-, tri)linear data on each cell, and most visualization programs can in fact only visualize this kind of data. That's good enough if you work with (bi-, tri)linear finite elements, in which case what you get to see is exactly what has been computed. On the other hand, if you work with (bi-, tri)quadratic elements, then what is written into the output file is just a (bi-, tri)linear interpolation onto the current mesh, i.e., only the values at the vertices. If this is not good enough, you can, for example, specify n_subdivisions equal to 2 to plot the solution on a once- refined mesh, or if set to 3, on a mesh where each cell is represented by 3-by-3 patches. On each of these smaller patches, given the limitations of output formats, the data is still linearly interpolated, but a linear interpolation of quadratic data on a finer mesh is still a better representation of the actual quadratic surface than on the original mesh. In other words, using this parameter can not help you plot the solution exactly, but it can get you closer if you use finite elements of higher polynomial degree.
Note
Specifying n_subdivisions>1 is useful when using higher order finite elements, but in general it does not actually result in the visualization showing higher order polynomial surfaces – rather, you just get a (bi-, tri-)linear interpolation of that higher order surface on a finer mesh. However, when outputting the solution in the VTK and VTU file formats via DataOutInterface::write_vtk() or DataOutInterface::write_vtu() (where DataOutInterface is a base class of the current class) as we often do in the tutorials, you can provide a set of flags via the DataOutBase::VtkFlags structure that includes the DataOutBase::VtkFlags::write_higher_order_cells flag. When set, the subdivisions produced by this function will be interpreted as support points for a higher order polynomial that will then actually be visualized as such. This is shown in step-11, for example. It is worth noting, however, that this requires a sufficiently new version of one of the VTK-based visualization programs.

Definition at line 1085 of file data_out.cc.

◆ build_patches() [2/3]

template<int dim, typename DoFHandlerType >
void DataOut< dim, DoFHandlerType >::build_patches ( const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  mapping,
const unsigned int  n_subdivisions = 0,
const CurvedCellRegion  curved_region = curved_boundary 
)
virtual

Same as above, except that the additional first parameter defines a mapping that is to be used in the generation of output. If n_subdivisions>1, the points interior of subdivided patches which originate from cells at the boundary of the domain can be computed using the mapping, i.e., a higher order mapping leads to a representation of a curved boundary by using more subdivisions. Some mappings like MappingQEulerian result in curved cells in the interior of the domain. The same is true if you have attached a manifold description to the cells of a triangulation (see Manifolds for more information). However, there is no easy way to query the mapping or manifold whether it really does lead to curved cells. Thus the last argument curved_region takes one of three values resulting in no curved cells at all, curved cells at the boundary (default) or curved cells in the whole domain. For more information about these three options, see the CurvedCellRegion enum's description.

Even for non-curved cells, the mapping argument can be used for Eulerian mappings (see class MappingQ1Eulerian) where a mapping is used not only to determine the position of points interior to a cell, but also of the vertices. It offers an opportunity to watch the solution on a deformed triangulation on which the computation was actually carried out, even if the mesh is internally stored in its undeformed configuration and the deformation is only tracked by an additional vector that holds the deformation of each vertex.

Definition at line 1101 of file data_out.cc.

◆ build_patches() [3/3]

template<int dim, typename DoFHandlerType >
void DataOut< dim, DoFHandlerType >::build_patches ( const hp::MappingCollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  mapping,
const unsigned int  n_subdivisions = 0,
const CurvedCellRegion  curved_region = curved_boundary 
)
virtual

Same as above, but for hp::MappingCollection.

Definition at line 1118 of file data_out.cc.

◆ set_cell_selection() [1/2]

template<int dim, typename DoFHandlerType >
void DataOut< dim, DoFHandlerType >::set_cell_selection ( const std::function< cell_iterator(const Triangulation< dim, spacedim > &)> &  first_cell,
const std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> &  next_cell 
)

A function that allows selecting for which cells output should be generated. This function takes two arguments, both std::function objects that can be used what the first cell on which output is generated is supposed to be, and what given one cell the next function is supposed to be. Through these function objects, it is possible to select a subset of cells on which output should be produced (e.g., only selecting those cells that belong to a part of the domain – say, the fluid domain in a code such as step-46), or to completely change where output is produced (e.g., to produce output on non-active cells of a multigrid hierarchy or if the finest level of a mesh is so fine that generating graphical output would lead to an overwhelming amount of data).

Parameters
[in]first_cellA function object that takes as argument the triangulation this class works on and that should return the first cell on which output should be generated.
[in]next_cellA function object that takes as arguments the triangulation as well as the last cell on which output was generated, and that should return the next cell on which output should be generated. If there is no next cell, i.e., if the input argument to the next_cell function object is the last cell on which output is to be generated, then next_cell must return triangulation.end().

These function objects are not difficult to write, but also not immediately obvious. As a consequence, there is a second variation of this function that takes a IteratorFilter argument and generates the corresponding functions itself.

Note
This function is also called in the constructor of this class, where the default behavior is set. By default, this class will select all locally owned and active cells for output.
If you have cell data (in contrast to nodal, or dof, data) such as error indicators, then you must make sure that the first_cell and next_cell function objects only walk over active cells, since cell data cannot be interpolated to a coarser cell. If you do have cell data and use this pair of functions and they return a non-active cell, then an exception will be thrown.

Definition at line 1309 of file data_out.cc.

◆ set_cell_selection() [2/2]

template<int dim, typename DoFHandlerType >
void DataOut< dim, DoFHandlerType >::set_cell_selection ( const FilteredIterator< cell_iterator > &  filtered_iterator)

A variation of the previous function that selects a subset of all cells for output based on the filter encoded in the FilteredIterator object given as argument. A typical way to generate the argument is via the make_filtered_iterator() function.

Alternatively, since FilteredIterator objects can be created from just a predicate (i.e., a function object that returns a bool), it is possible to call this function with just a lambda function, which will then automatically be converted to a FilteredIterator object. For example, the following piece of code works:

DataOut<dim> data_out;
[](const typename Triangulation<dim>::cell_iterator &cell) {
return (cell->is_active() && cell->subdomain_id() == 0);
});
void set_cell_selection(const std::function< cell_iterator(const Triangulation< dim, spacedim > &)> &first_cell, const std::function< cell_iterator(const Triangulation< dim, spacedim > &, const cell_iterator &)> &next_cell)
Definition: data_out.cc:1309

In this case, the lambda function selects all of those cells that are active and whose subdomain id is zero. These will then be the only cells on which output is generated.

Note
Not all filters will result in subsets of cells for which output can actually be generated. For example, if you are working on parallel meshes where data is only available on some cells, then you better make sure that your filtered_iterator only loops over the locally owned cells; likewise, in most cases you will probably only want to work on active cells since this is where the solution actually lives. In particular, if you have added vectors that represent data defined on cells (instead of nodal data), then you can not generate output on non-active cells and your iterator filter should reflect this.

Definition at line 1323 of file data_out.cc.

◆ get_cell_selection()

template<int dim, typename DoFHandlerType >
const std::pair< typename DataOut< dim, DoFHandlerType >::FirstCellFunctionType, typename DataOut< dim, DoFHandlerType >::NextCellFunctionType > DataOut< dim, DoFHandlerType >::get_cell_selection

Return the two function objects that are in use for determining the first and the next cell as set by set_cell_selection().

Definition at line 1365 of file data_out.cc.

◆ first_cell()

template<int dim, typename DoFHandlerType >
DataOut< dim, DoFHandlerType >::cell_iterator DataOut< dim, DoFHandlerType >::first_cell
virtual

Return the first cell which we want output for. The default implementation returns the first active cell, but you might want to return other cells in a derived class.

Deprecated:
Use the set_cell_selection() function instead.

Definition at line 1374 of file data_out.cc.

◆ next_cell()

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
DataOut< dim, DoFHandlerType >::cell_iterator DataOut< dim, DoFHandlerType >::next_cell ( const cell_iterator cell)
virtual

Return the next cell after cell which we want output for. If there are no more cells, any implementation of this function should return dof_handler->end().

The default implementation returns the next active cell, but you might want to return other cells in a derived class. Note that the default implementation assumes that the given cell is active, which is guaranteed as long as first_cell() is also used from the default implementation. Overloading only one of the two functions might not be a good idea.

Deprecated:
Use the set_cell_selection() function instead.

Definition at line 1383 of file data_out.cc.

◆ first_locally_owned_cell()

template<int dim, typename DoFHandlerType >
DataOut< dim, DoFHandlerType >::cell_iterator DataOut< dim, DoFHandlerType >::first_locally_owned_cell
privatevirtual

Return the first cell produced by the first_cell()/next_cell() function pair that is locally owned. If this object operates on a non-distributed triangulation, the result equals what first_cell() returns.

Deprecated:
Use the set_cell_selection() function instead.

Definition at line 1399 of file data_out.cc.

◆ next_locally_owned_cell()

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
DataOut< dim, DoFHandlerType >::cell_iterator DataOut< dim, DoFHandlerType >::next_locally_owned_cell ( const cell_iterator cell)
privatevirtual

Return the next cell produced by the next_cell() function that is locally owned. If this object operates on a non-distributed triangulation, the result equals what first_cell() returns.

Deprecated:
Use the set_cell_selection() function instead.

Definition at line 1416 of file data_out.cc.

◆ build_one_patch()

template<int dim, typename DoFHandlerType >
void DataOut< dim, DoFHandlerType >::build_one_patch ( const std::pair< cell_iterator, unsigned int > *  cell_and_index,
internal::DataOutImplementation::ParallelData< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  scratch_data,
const unsigned int  n_subdivisions,
const CurvedCellRegion  curved_cell_region 
)
private

Build one patch. This function is called in a WorkStream context.

The first argument here is the iterator, the second the scratch data object. All following are tied to particular values when calling WorkStream::run(). The function does not take a CopyData object but rather allocates one on its own stack for memory access efficiency reasons.

Definition at line 88 of file data_out.cc.

◆ attach_dof_handler()

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::attach_dof_handler ( const DoFHandler< dim > &  )
inherited

Designate a dof handler to be used to extract geometry data and the mapping between nodes and node values. This call is not necessary if all added data vectors are supplemented with a DoFHandler argument.

This call is optional: If you add data vectors with specified DoFHandler object, then that contains all information needed to generate the output.

◆ attach_triangulation()

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::attach_triangulation ( const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  )
inherited

Designate a triangulation to be used to extract geometry data and the mapping between nodes and node values.

This call is optional: If you add data vectors with specified DoFHandler object, then that contains all information needed to generate the output. This call is useful when you only output cell vectors and no DoFHandler at all, in which case it provides the geometry.

◆ add_data_vector() [1/6]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector ( const VectorType &  data,
const std::vector< std::string > &  names,
const DataVectorType  type = type_automatic,
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretation = std::vector< DataComponentInterpretation::DataComponentInterpretation>() 
)
inherited

Add a data vector together with its name.

A pointer to the vector is stored, so you have to make sure the vector exists at that address at least as long as you call the write_* functions.

It is assumed that the vector has the same number of components as there are degrees of freedom in the dof handler, in which case it is assumed to be a vector storing nodal data; or the size may be the number of active cells on the present grid, in which case it is assumed to be a cell data vector. As the number of degrees of freedom and of cells is usually not equal, the function can determine itself which type of vector it is given. However, there are corner cases where this automatic determination does not work. One example is if you compute with piecewise constant elements and have a scalar solution, then there are as many cells as there are degrees of freedom (though they may be numbered differently). Another possibility is if you have a 1d mesh embedded in 2d space and the mesh consists of a closed curve of cells; in this case, there are as many nodes as there are cells, and when using a Q1 element you will have as many degrees of freedom as there are cells. In these cases, you can change the last argument of the function from its default value type_automatic to either type_dof_data or type_cell_data, depending on what the vector represents. Apart from such corner cases, you can leave the argument at its default value and let the function determine the type of the vector itself.

If it is a vector holding DoF data, the names given shall be one for each component of the underlying finite element. If it is a finite element composed of only one subelement, then there is another function following which takes a single name instead of a vector of names.

The data_component_interpretation argument contains information about how the individual components of output files that consist of more than one data set are to be interpreted.

For example, if one has a finite element for the Stokes equations in 2d, representing components (u,v,p), one would like to indicate that the first two, u and v, represent a logical vector so that later on when we generate graphical output we can hand them off to a visualization program that will automatically know to render them as a vector field, rather than as two separate and independent scalar fields.

The default value of this argument (i.e. an empty vector) corresponds is equivalent to a vector of values DataComponentInterpretation::component_is_scalar, indicating that all output components are independent scalar fields. However, if the given data vector represents logical vectors, you may pass a vector that contains values DataComponentInterpretation::component_is_part_of_vector. In the example above, one would pass in a vector with components (DataComponentInterpretation::component_is_part_of_vector, DataComponentInterpretation::component_is_part_of_vector, DataComponentInterpretation::component_is_scalar) for (u,v,p).

The names of a data vector shall only contain characters which are letters, underscore and a few other ones. Refer to the ExcInvalidCharacter exception declared in this class to see which characters are valid and which are not.

Note
The actual type for the vector argument may be any vector type from which FEValues can extract values on a cell using the FEValuesBase::get_function_values() function.
When working in parallel, the vector to be written needs to be ghosted with read access to all degrees of freedom on the locally owned cells, see the step-40 or step-37 tutorial programs for details, i.e., it might be necessary to call data.update_ghost_values().

Definition at line 739 of file data_out_dof_data.h.

◆ add_data_vector() [2/6]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector ( const VectorType &  data,
const std::string &  name,
const DataVectorType  type = type_automatic,
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretation = std::vector< DataComponentInterpretation::DataComponentInterpretation>() 
)
inherited

This function is an abbreviation to the above one (see there for a discussion of the various arguments), intended for use with finite elements that are not composed of subelements. In this case, only one name per data vector needs to be given, which is what this function takes. It simply relays its arguments after a conversion of the name to a vector of strings, to the other add_data_vector() function above.

If data is a vector with multiple components this function will generate distinct names for all components by appending an underscore and the number of each component to name

The actual type for the template argument may be any vector type from which FEValues can extract values on a cell using the FEValuesBase::get_function_values() function.

Definition at line 765 of file data_out_dof_data.h.

◆ add_data_vector() [3/6]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector ( const DoFHandler< dim > &  dof_handler,
const VectorType &  data,
const std::vector< std::string > &  names,
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretation = std::vector< DataComponentInterpretation::DataComponentInterpretation>() 
)
inherited

This function is an extension of the above one (see there for a discussion of the arguments except the first one) and allows to set a vector with its own DoFHandler object. This DoFHandler needs to be compatible with the other DoFHandler objects assigned with calls to add_data_vector or attach_dof_handler, in the sense that all of the DoFHandler objects need to be based on the same triangulation. This function allows you to export data from multiple DoFHandler objects that describe different solution components. An example of using this function is given in step-61.

Since this function takes a DoFHandler object and hence naturally represents dof data, the data vector type argument present in the other methods above is not necessary.

Definition at line 790 of file data_out_dof_data.h.

◆ add_data_vector() [4/6]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector ( const DoFHandler< dim > &  dof_handler,
const VectorType &  data,
const std::string &  name,
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretation = std::vector< DataComponentInterpretation::DataComponentInterpretation>() 
)
inherited

This function is an abbreviation of the function above with only a scalar dof_handler given and a single data name.

Definition at line 805 of file data_out_dof_data.h.

◆ add_data_vector() [5/6]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector ( const VectorType &  data,
const DataPostprocessor< DoFHandlerType::space_dimension > &  data_postprocessor 
)
inherited

This function is an alternative to the above ones, allowing the output of derived quantities instead of the given data. This conversion has to be done in a class derived from DataPostprocessor. This function is used in step-29. Other uses are shown in step-32 and step-33.

The names for these derived quantities are provided by the data_postprocessor argument. Likewise, the data_component_interpretation argument of the other add_data_vector() functions is provided by the data_postprocessor argument. As only data of type type_dof_data can be transformed, this type is also known implicitly and does not have to be given.

Note
The actual type for the vector argument may be any vector type from which FEValues can extract values on a cell using the FEValuesBase::get_function_values() function.
The DataPostprocessor object (i.e., in reality the object of your derived class) has to live until the DataOut object is destroyed as the latter keeps a pointer to the former and will complain if the object pointed to is destroyed while the latter still has a pointer to it. If both the data postprocessor and DataOut objects are local variables of a function (as they are, for example, in step-29), then you can avoid this error by declaring the data postprocessor variable before the DataOut variable as objects are destroyed in reverse order of declaration.

Definition at line 841 of file data_out_dof_data.h.

◆ add_data_vector() [6/6]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector ( const DoFHandler< dim > &  dof_handler,
const VectorType &  data,
const DataPostprocessor< DoFHandlerType::space_dimension > &  data_postprocessor 
)
inherited

Same function as above, but with a DoFHandler object that does not need to coincide with the DoFHandler initially set. Note that the postprocessor can only read data from the given DoFHandler and solution vector, not other solution vectors or DoFHandlers.

◆ add_mg_data_vector() [1/2]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_mg_data_vector ( const DoFHandler< dim > &  dof_handler,
const MGLevelObject< VectorType > &  data,
const std::vector< std::string > &  names,
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretation = std::vector< DataComponentInterpretation::DataComponentInterpretation>() 
)
inherited

Add a multilevel data vector.

This function adds the vector-valued multilevel vector data in the form of a vector on each level that belongs to the DoFHandler dof_handler to the graphical output. This function is typically used in conjunction with a call to set_cell_selection() that selects cells on a specific level and not the active cells (the default).

A vector data can be obtained in several ways, for example by using Multigrid::solution or Multigrid::defect during or after a multigrid cycle or by interpolating a solution via MGTransferMatrixFree::interpolate_to_mg().

The handling of names and data_component_interpretation is identical to the add_data_vector() function.

◆ add_mg_data_vector() [2/2]

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_mg_data_vector ( const DoFHandler< dim > &  dof_handler,
const MGLevelObject< VectorType > &  data,
const std::string &  name 
)
inherited

Scalar version of the function above.

◆ clear_data_vectors()

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::clear_data_vectors ( )
inherited

Release the pointers to the data vectors. This allows output of a new set of vectors without supplying the DoF handler again. Therefore, the DataOut object can be used in an algebraic context. Note that besides the data vectors also the patches already computed are deleted.

◆ clear_input_data_references()

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::clear_input_data_references ( )
inherited

Release pointers to all input data elements, i.e. pointers to data vectors and to the DoF handler object. This function may be useful when you have called the build_patches function of derived class, since then the patches are built and the input data is no more needed, nor is there a need to reference it. You can then output the patches detached from the main thread and need not make sure anymore that the DoF handler object and vectors must not be deleted before the output thread is finished.

◆ merge_patches()

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::merge_patches ( const DataOut_DoFData< DoFHandlerType2, patch_dim, patch_space_dim > &  source,
const Point< patch_space_dim > &  shift = Point<patch_space_dim>() 
)
inherited

This function can be used to merge the patches that were created using the build_patches function of the object given as argument into the list of patches created by this object. This is sometimes handy if one has, for example, a domain decomposition algorithm where each block is represented by a DoFHandler of its own, but one wants to output the solution on all the blocks at the same time.

For this to work, the given argument and this object need to have the same number of output vectors, and they need to use the same number of subdivisions per patch. The output will probably look rather funny if patches in both objects overlap in space.

If you call build_patches() for this object after merging in patches, the previous state is overwritten, and the merged-in patches are lost.

The second parameter allows to shift each node of the patches in the object passed in in the first parameter by a certain amount. This is sometimes useful to generate "exploded" views of a collection of blocks.

This function will fail if either this or the other object did not yet set up any patches.

Definition at line 941 of file data_out_dof_data.h.

◆ clear()

virtual void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::clear ( )
virtualinherited

Release the pointers to the data vectors and the DoF handler. You have to set all data entries again using the add_data_vector() function. The pointer to the dof handler is cleared as well, along with all other data. In effect, this function resets everything to a virgin state.

◆ memory_consumption()

std::size_t DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::memory_consumption ( ) const
inherited

Determine an estimate for the memory consumption (in bytes) of this object.

◆ get_patches()

virtual const std::vector< Patch > & DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::get_patches ( ) const
overrideprotectedvirtualinherited

Function by which the base class's functions get to know what patches they shall write to a file.

Implements DataOutInterface< patch_dim, patch_space_dim >.

◆ get_dataset_names()

virtual std::vector< std::string > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::get_dataset_names ( ) const
overrideprotectedvirtualinherited

Virtual function through which the names of data sets are obtained by the output functions of the base class.

Implements DataOutInterface< patch_dim, patch_space_dim >.

◆ get_fes()

std::vector< std::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::get_fes ( ) const
protectedinherited

Extracts the finite elements stored in the dof_data object, including a dummy object of FE_DGQ<dim>(0) in case only the triangulation is used.

◆ get_nonscalar_data_ranges()

virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::get_nonscalar_data_ranges ( ) const
overrideprotectedvirtualinherited

Overload of the respective DataOutInterface::get_nonscalar_data_ranges() function. See there for a more extensive documentation.

Reimplemented from DataOutInterface< patch_dim, patch_space_dim >.

◆ add_data_vector_internal()

void DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::add_data_vector_internal ( const DoFHandler< dim > *  dof_handler,
const VectorType &  data,
const std::vector< std::string > &  names,
const DataVectorType  type,
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretation,
const bool  deduce_output_names 
)
privateinherited

Common function called by the four public add_data_vector methods.

◆ write_dx()

void DataOutInterface< dim, spacedim >::write_dx ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in OpenDX format. See DataOutBase::write_dx.

Definition at line 2576 of file data_out_base.cc.

◆ write_eps()

void DataOutInterface< dim, spacedim >::write_eps ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in EPS format. See DataOutBase::write_eps.

Definition at line 2583 of file data_out_base.cc.

◆ write_gmv()

void DataOutInterface< dim, spacedim >::write_gmv ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in GMV format. See DataOutBase::write_gmv.

Definition at line 2590 of file data_out_base.cc.

◆ write_gnuplot()

void DataOutInterface< dim, spacedim >::write_gnuplot ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in GNUPLOT format. See DataOutBase::write_gnuplot.

Definition at line 2597 of file data_out_base.cc.

◆ write_povray()

void DataOutInterface< dim, spacedim >::write_povray ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in POVRAY format. See DataOutBase::write_povray.

Definition at line 2604 of file data_out_base.cc.

◆ write_tecplot()

void DataOutInterface< dim, spacedim >::write_tecplot ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in Tecplot format. See DataOutBase::write_tecplot.

Definition at line 2611 of file data_out_base.cc.

◆ write_ucd()

void DataOutInterface< dim, spacedim >::write_ucd ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in UCD format for AVS. See DataOutBase::write_ucd.

Definition at line 2618 of file data_out_base.cc.

◆ write_vtk()

void DataOutInterface< dim, spacedim >::write_vtk ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in Vtk format. See DataOutBase::write_vtk.

Note
VTK is a legacy format and has largely been supplanted by the VTU format (an XML-structured version of VTK). In particular, VTU allows for the compression of data and consequently leads to much smaller file sizes that equivalent VTK files for large files. Since all visualization programs that support VTK also support VTU, you should consider using the latter file format instead, by using the write_vtu() function.

Definition at line 2632 of file data_out_base.cc.

◆ write_vtu()

void DataOutInterface< dim, spacedim >::write_vtu ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in Vtu (VTK's XML) format. See DataOutBase::write_vtu.

Some visualization programs, such as ParaView, can read several separate VTU files to parallelize visualization. In that case, you need a .pvtu file that describes which VTU files form a group. The DataOutInterface::write_pvtu_record() function can generate such a centralized record. Likewise, DataOutInterface::write_visit_record() does the same for older versions of VisIt (although VisIt can also read pvtu records since version 2.5.1). Finally, DataOutInterface::write_pvd_record() can be used to group together the files that jointly make up a time dependent simulation.

Definition at line 2649 of file data_out_base.cc.

◆ write_vtu_in_parallel()

void DataOutInterface< dim, spacedim >::write_vtu_in_parallel ( const std::string &  filename,
const MPI_Comm comm 
) const
inherited

Collective MPI call to write the solution from all participating nodes (those in the given communicator) to a single compressed .vtu file on a shared file system. The communicator can be a sub communicator of the one used by the computation. This routine uses MPI I/O to achieve high performance on parallel filesystems. Also see DataOutInterface::write_vtu().

Definition at line 2660 of file data_out_base.cc.

◆ write_pvtu_record()

void DataOutInterface< dim, spacedim >::write_pvtu_record ( std::ostream &  out,
const std::vector< std::string > &  piece_names 
) const
inherited

Some visualization programs, such as ParaView, can read several separate VTU files that all form part of the same simulation, in order to parallelize visualization. In that case, you need a .pvtu file that describes which VTU files (written, for example, through the DataOutInterface::write_vtu() function) form a group. The current function can generate such a centralized record.

The central record file generated by this function contains a list of (scalar or vector) fields that describes which fields can actually be found in the individual files that comprise the set of parallel VTU files along with the names of these files. This function gets the names and types of fields through the get_dataset_names() and get_nonscalar_data_ranges() functions of this class. The second argument to this function specifies the names of the files that form the parallel set.

Note
Use DataOutBase::write_vtu() and DataOutInterface::write_vtu() for writing each piece. Also note that only one parallel process needs to call the current function, listing the names of the files written by all parallel processes.
The use of this function is explained in step-40.
In order to tell Paraview to group together multiple pvtu files that each describe one time step of a time dependent simulation, see the DataOutBase::write_pvd_record() function.
Older versions of VisIt (before 2.5.1), can not read pvtu records. However, it can read visit records as written by the write_visit_record() function.

Definition at line 2697 of file data_out_base.cc.

◆ write_vtu_with_pvtu_record()

std::string DataOutInterface< dim, spacedim >::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
inherited

This function writes several .vtu files and a .pvtu record in parallel and constructs the filenames automatically. It is a combination of DataOutInterface::write_vtu() or DataOutInterface::write_vtu_in_parallel(), and DataOutInterface::write_pvtu_record().

For example, running write_vtu_with_pvtu_record("output/", "solution", 3, comm, 4, 2) on 10 processes generates the files

output/solution_0003.0.vtu
output/solution_0003.1.vtu
output/solution_0003.pvtu

where the .0.vtu file contains the output of the first half of the processes grouped together, and the .1.vtu the data from the remaining half.

A specified directory and a filename_without_extension form the first part of the filename. The filename is then extended with a counter labeling the current timestep/iteration/etc., the processor ID, and finally the .vtu/.pvtu ending. Since the number of timesteps to be written depends on the application, the number of digits to be reserved in the filename can be specified as parameter n_digits_for_counter, and the number is not padded with leading zeros if this parameter is left at its default value numbers::invalid_unsigned_int. If more than one file identifier is needed (e.g. time step number and iteration counter of solver), the last identifier is used as counter, while all other identifiers have to be added to filename_without_extension when calling this function.

In a parallel setting, several files are typically written per time step. The number of files written in parallel depends on the number of MPI processes (see parameter mpi_communicator), and a specified number of n_groups with default value 0. The background is that VTU file output supports grouping files from several CPUs into a given number of files using MPI I/O when writing on a parallel filesystem. The default value of n_groups is 0, meaning that every MPI rank will write one file. A value of 1 will generate one big file containing the solution over the whole domain, while a larger value will create n_groups files (but not more than there are MPI ranks).

Note that only one processor needs to generate the .pvtu file, where processor zero is chosen to take over this job.

The return value is the filename of the centralized file for the pvtu record.

Note
The code simply combines the strings directory and filename_without_extension, i.e., the user has to make sure that directory contains a trailing character, e.g. "/", that separates the directory from the filename.
Use an empty string "" for the first argument if output is to be written in the current working directory.

Definition at line 2759 of file data_out_base.cc.

◆ write_svg()

void DataOutInterface< dim, spacedim >::write_svg ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in SVG format. See DataOutBase::write_svg.

Definition at line 2772 of file data_out_base.cc.

◆ write_deal_II_intermediate()

void DataOutInterface< dim, spacedim >::write_deal_II_intermediate ( std::ostream &  out) const
inherited

Obtain data through get_patches() and write it to out in deal.II intermediate format. See DataOutBase::write_deal_II_intermediate.

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.

Definition at line 2785 of file data_out_base.cc.

◆ create_xdmf_entry() [1/2]

XDMFEntry DataOutInterface< dim, spacedim >::create_xdmf_entry ( const DataOutBase::DataOutFilter data_filter,
const std::string &  h5_filename,
const double  cur_time,
const MPI_Comm comm 
) const
inherited

Create an XDMFEntry based on the data in the data_filter. This assumes the mesh and solution data were written to a single file. See write_xdmf_file() for an example of usage.

Definition at line 2793 of file data_out_base.cc.

◆ create_xdmf_entry() [2/2]

XDMFEntry DataOutInterface< dim, spacedim >::create_xdmf_entry ( const DataOutBase::DataOutFilter data_filter,
const std::string &  h5_mesh_filename,
const std::string &  h5_solution_filename,
const double  cur_time,
const MPI_Comm comm 
) const
inherited

Create an XDMFEntry based on the data in the data_filter. This assumes the mesh and solution data were written to separate files. See write_xdmf_file() for an example of usage.

Definition at line 2804 of file data_out_base.cc.

◆ write_xdmf_file()

void DataOutInterface< dim, spacedim >::write_xdmf_file ( const std::vector< XDMFEntry > &  entries,
const std::string &  filename,
const MPI_Comm comm 
) const
inherited

Write an XDMF file based on the provided vector of XDMFEntry objects. Below is an example of how to use this function with HDF5 and the DataOutFilter:

DataOutBase::DataOutFilter data_filter(flags);
std::vector<XDMFEntry> xdmf_entries;
// Filter the data and store it in data_filter
data_out.write_filtered_data(data_filter);
// Write the filtered data to HDF5
data_out.write_hdf5_parallel(data_filter, "solution.h5", MPI_COMM_WORLD);
// Create an XDMF entry detailing the HDF5 file
auto new_xdmf_entry = data_out.create_xdmf_entry(data_filter,
"solution.h5",
simulation_time,
MPI_COMM_WORLD);
// Add the XDMF entry to the list
xdmf_entries.push_back(new_xdmf_entry);
// Create an XDMF file from all stored entries
data_out.write_xdmf_file(xdmf_entries, "solution.xdmf", MPI_COMM_WORLD);

Definition at line 2835 of file data_out_base.cc.

◆ write_hdf5_parallel() [1/2]

void DataOutInterface< dim, spacedim >::write_hdf5_parallel ( const DataOutBase::DataOutFilter data_filter,
const std::string &  filename,
const MPI_Comm comm 
) const
inherited

Write the data in data_filter to a single HDF5 file containing both the mesh and solution values. Below is an example of how to use this function with the DataOutFilter:

DataOutBase::DataOutFilter data_filter(flags);
// Filter the data and store it in data_filter
data_out.write_filtered_data(data_filter);
// Write the filtered data to HDF5
data_out.write_hdf5_parallel(data_filter, "solution.h5", MPI_COMM_WORLD);

Definition at line 2854 of file data_out_base.cc.

◆ write_hdf5_parallel() [2/2]

void DataOutInterface< dim, spacedim >::write_hdf5_parallel ( const DataOutBase::DataOutFilter data_filter,
const bool  write_mesh_file,
const std::string &  mesh_filename,
const std::string &  solution_filename,
const MPI_Comm comm 
) const
inherited

Write the data in data_filter to HDF5 file(s). If write_mesh_file is false, the mesh data will not be written and the solution file will contain only the solution values. If write_mesh_file is true and the filenames are the same, the resulting file will contain both mesh data and solution values.

Definition at line 2866 of file data_out_base.cc.

◆ write_filtered_data()

void DataOutInterface< dim, spacedim >::write_filtered_data ( DataOutBase::DataOutFilter filtered_data) const
inherited

DataOutFilter is an intermediate data format that reduces the amount of data that will be written to files. The object filled by this function can then later be used again to write data in a concrete file format; see, for example, DataOutBase::write_hdf5_parallel().

Definition at line 2879 of file data_out_base.cc.

◆ write()

void DataOutInterface< dim, spacedim >::write ( std::ostream &  out,
const DataOutBase::OutputFormat  output_format = DataOutBase::default_format 
) const
inherited

Write data and grid to out according to the given data format. This function simply calls the appropriate write_* function. If no output format is requested, the default_format is written.

An error occurs if no format is provided and the default format is default_format.

Definition at line 2891 of file data_out_base.cc.

◆ set_default_format()

void DataOutInterface< dim, spacedim >::set_default_format ( const DataOutBase::OutputFormat  default_format)
inherited

Set the default format. The value set here is used anytime, output for format default_format is requested.

Definition at line 2900 of file data_out_base.cc.

◆ set_flags()

void DataOutInterface< dim, spacedim >::set_flags ( const FlagType &  flags)
inherited

Set the flags to be used for output. This method expects flags to be a member of one of the child classes of OutputFlagsBase.

Definition at line 2909 of file data_out_base.cc.

◆ default_suffix()

std::string DataOutInterface< dim, spacedim >::default_suffix ( const DataOutBase::OutputFormat  output_format = DataOutBase::default_format) const
inherited

A function that returns the same string as the respective function in the base class does; the only exception being that if the parameter is omitted, then the value for the present default format is returned, i.e. the correct suffix for the format that was set through set_default_format() or parse_parameters() before calling this function.

Definition at line 2920 of file data_out_base.cc.

◆ declare_parameters()

void DataOutInterface< dim, spacedim >::declare_parameters ( ParameterHandler prm)
staticinherited

Declare parameters for all output formats by declaring subsections within the parameter file for each output format and call the respective declare_parameters functions of the flag classes for each output format.

Some of the declared subsections may not contain entries, if the respective format does not export any flags.

Note that the top-level parameters denoting the number of subdivisions per patch and the output format are not declared, since they are only passed to virtual functions and are not stored inside objects of this type. You have to declare them yourself.

Definition at line 2938 of file data_out_base.cc.

◆ parse_parameters()

void DataOutInterface< dim, spacedim >::parse_parameters ( ParameterHandler prm)
inherited

Read the parameters declared in declare_parameters() and set the flags for the output formats accordingly.

The flags thus obtained overwrite all previous contents of the flag objects as default-constructed or set by the set_flags() function.

Definition at line 2948 of file data_out_base.cc.

◆ validate_dataset_names()

void DataOutInterface< dim, spacedim >::validate_dataset_names
protectedinherited

Validate that the names of the datasets returned by get_dataset_names() and get_nonscalar_data_ranges() are valid. This currently consists of checking that names are not used more than once. If an invalid state is encountered, an Assert() will be triggered in debug mode.

Definition at line 3008 of file data_out_base.cc.

Member Data Documentation

◆ spacedim

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
constexpr unsigned int DataOut< dim, DoFHandlerType >::spacedim = DoFHandlerType::space_dimension
staticconstexpr

Definition at line 156 of file data_out.h.

◆ first_cell_function

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
std::function<cell_iterator(const Triangulation<dim, spacedim> &)> DataOut< dim, DoFHandlerType >::first_cell_function
private

A function object that is used to select what the first cell is going to be on which to generate graphical output. See the set_cell_selection() function for more information.

Definition at line 468 of file data_out.h.

◆ next_cell_function

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
std::function<cell_iterator(const Triangulation<dim, spacedim> &, const cell_iterator &)> DataOut< dim, DoFHandlerType >::next_cell_function
private

A function object that is used to select what the next cell is going to be on which to generate graphical output, given a previous cell. See the set_cell_selection() function for more information.

Definition at line 477 of file data_out.h.

◆ triangulation

SmartPointer<const Triangulation<DoFHandlerType::dimension, DoFHandlerType::space_dimension> > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::triangulation
protectedinherited

Pointer to the triangulation object.

Definition at line 972 of file data_out_dof_data.h.

◆ dofs

SmartPointer<const DoFHandler< dim > > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::dofs
protectedinherited

Pointer to the optional handler object.

Definition at line 977 of file data_out_dof_data.h.

◆ dof_data

std::vector<std::shared_ptr<internal::DataOutImplementation::DataEntryBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension> > > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::dof_data
protectedinherited

List of data elements with vectors of values for each degree of freedom.

Definition at line 985 of file data_out_dof_data.h.

◆ cell_data

std::vector<std::shared_ptr<internal::DataOutImplementation::DataEntryBase< DoFHandlerType::dimension, DoFHandlerType::space_dimension> > > DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::cell_data
protectedinherited

List of data elements with vectors of values for each cell.

Definition at line 993 of file data_out_dof_data.h.

◆ patches

std::vector<Patch> DataOut_DoFData< DoFHandler< dim > , patch_dim, patch_space_dim >::patches
protectedinherited

This is a list of patches that is created each time build_patches() is called. These patches are used in the output routines of the base classes.

Definition at line 1000 of file data_out_dof_data.h.

◆ default_subdivisions

unsigned int DataOutInterface< dim, spacedim >::default_subdivisions
protectedinherited

The default number of subdivisions for patches. This is filled by parse_parameters() and should be obeyed by build_patches() in derived classes.

Definition at line 3016 of file data_out_base.h.

◆ default_fmt

DataOutBase::OutputFormat DataOutInterface< dim, spacedim >::default_fmt
privateinherited

Standard output format. Use this format, if output format default_format is requested. It can be changed by the set_format function or in a parameter file.

Definition at line 3024 of file data_out_base.h.

◆ dx_flags

DataOutBase::DXFlags DataOutInterface< dim, spacedim >::dx_flags
privateinherited

Flags to be used upon output of OpenDX data. Can be changed by using the set_flags function.

Definition at line 3030 of file data_out_base.h.

◆ ucd_flags

DataOutBase::UcdFlags DataOutInterface< dim, spacedim >::ucd_flags
privateinherited

Flags to be used upon output of UCD data. Can be changed by using the set_flags function.

Definition at line 3036 of file data_out_base.h.

◆ gnuplot_flags

DataOutBase::GnuplotFlags DataOutInterface< dim, spacedim >::gnuplot_flags
privateinherited

Flags to be used upon output of GNUPLOT data. Can be changed by using the set_flags function.

Definition at line 3042 of file data_out_base.h.

◆ povray_flags

DataOutBase::PovrayFlags DataOutInterface< dim, spacedim >::povray_flags
privateinherited

Flags to be used upon output of POVRAY data. Can be changed by using the set_flags function.

Definition at line 3048 of file data_out_base.h.

◆ eps_flags

DataOutBase::EpsFlags DataOutInterface< dim, spacedim >::eps_flags
privateinherited

Flags to be used upon output of EPS data in one space dimension. Can be changed by using the set_flags function.

Definition at line 3054 of file data_out_base.h.

◆ gmv_flags

DataOutBase::GmvFlags DataOutInterface< dim, spacedim >::gmv_flags
privateinherited

Flags to be used upon output of gmv data in one space dimension. Can be changed by using the set_flags function.

Definition at line 3060 of file data_out_base.h.

◆ tecplot_flags

DataOutBase::TecplotFlags DataOutInterface< dim, spacedim >::tecplot_flags
privateinherited

Flags to be used upon output of Tecplot data in one space dimension. Can be changed by using the set_flags function.

Definition at line 3066 of file data_out_base.h.

◆ vtk_flags

DataOutBase::VtkFlags DataOutInterface< dim, spacedim >::vtk_flags
privateinherited

Flags to be used upon output of vtk data in one space dimension. Can be changed by using the set_flags function.

Definition at line 3072 of file data_out_base.h.

◆ svg_flags

DataOutBase::SvgFlags DataOutInterface< dim, spacedim >::svg_flags
privateinherited

Flags to be used upon output of svg data in one space dimension. Can be changed by using the set_flags function.

Definition at line 3078 of file data_out_base.h.

◆ deal_II_intermediate_flags

DataOutBase::Deal_II_IntermediateFlags DataOutInterface< dim, spacedim >::deal_II_intermediate_flags
privateinherited

Flags to be used upon output of deal.II intermediate data in one space dimension. Can be changed by using the set_flags function.

Definition at line 3084 of file data_out_base.h.


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