Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | Public Member Functions | Private Member Functions | 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
 
- Public Types inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension >
enum  DataVectorType
 
using cell_iterator = typename Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension >::cell_iterator
 

Public Member Functions

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 cell_iterator first_cell ()
 
virtual cell_iterator next_cell (const cell_iterator &cell)
 
- Public Member Functions inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension >
 DataOut_DoFData ()
 
virtual ~DataOut_DoFData () override
 
void attach_dof_handler (const DoFHandlerType &)
 
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 DoFHandlerType &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 DoFHandlerType &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 DoFHandlerType &dof_handler, const VectorType &data, const DataPostprocessor< DoFHandlerType::space_dimension > &data_postprocessor)
 
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
 
- Public Member Functions inherited from DataOutInterface< patch_dim, patch_space_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
 
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
 

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)
 

Additional Inherited Members

- Static Public Member Functions inherited from DataOutInterface< patch_dim, patch_space_dim >
static void declare_parameters (ParameterHandler &prm)
 
- Protected Types inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension >
using Patch = ::DataOutBase::Patch< patch_dim, patch_space_dim >
 
- Protected Member Functions inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension >
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
 
- Protected Member Functions inherited from DataOutInterface< patch_dim, patch_space_dim >
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string > > get_vector_data_ranges () const
 
void validate_dataset_names () const
 
- Protected Attributes inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension >
SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > triangulation
 
SmartPointer< const DoFHandlerType > dofs
 
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType > > > dof_data
 
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType > > > cell_data
 
std::vector< Patchpatches
 
- Protected Attributes inherited from DataOutInterface< patch_dim, patch_space_dim >
unsigned int default_subdivisions
 

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. 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 three, and so on. 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 they 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 in parallel programs such as the step-18 example program), 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 functions first_cell() and next_cell(). By default, they return the first active cell, and the next active cell, respectively. Since they are virtual functions, you can write your own class derived from DataOut in which you overload these two functions to select other cells for output. 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 first_cell() and next_cell() 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). Once you derive your own class, you would just create an object of this type instead of an object of type DataOut, and everything else will remain the same.

The two functions are not constant, so you may store information within your derived class about the last accessed cell. This is useful if the information of the last cell which was accessed is not sufficient to determine the next one.

There is one caveat, however: if you have cell data (in contrast to nodal, or dof, data) such as error indicators, then you must make sure that first_cell() and next_cell() 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.

Precondition
This class only makes sense if the first template argument, dim equals the dimension of the DoFHandler type given as the second template argument, i.e., if dim == DoFHandlerType::dimension. This redundancy is a historical relic from the time where the library had only a single DoFHandler class and this class consequently only a single template argument.
Author
Wolfgang Bangerth, 1999

Definition at line 160 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 172 of file data_out.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 faces).

Definition at line 191 of file data_out.h.

Member Function Documentation

◆ build_patches() [1/2]

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 point for a higher order polynomial that will then actually be visualized as such. On the other hand, this requires a sufficiently new version of one of the VTK-based visualization programs.

Definition at line 423 of file data_out.cc.

◆ build_patches() [2/2]

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.

Todo:
The mapping argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler.

Definition at line 435 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.

Definition at line 610 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, dofs->end() shall be returned.

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.

Definition at line 619 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.

Definition at line 635 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.

Definition at line 652 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 70 of file data_out.cc.


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