Reference documentation for deal.II version 9.2.0
|
#include <deal.II/numerics/data_out.h>
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 &)> |
Public Types inherited from DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension > | |
enum | DataVectorType |
using | cell_iterator = typename Triangulation< DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension >::cell_iterator |
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) |
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, NextCellFunctionType > | get_cell_selection () const |
virtual cell_iterator | first_cell () |
virtual cell_iterator | next_cell (const cell_iterator &cell) |
Public Member Functions inherited from DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension > | |
DataOut_DoFData () | |
virtual | ~DataOut_DoFData () override |
void | attach_dof_handler (const DoFHandler< dim > &) |
void | attach_triangulation (const Triangulation< DoFHandler< dim > ::dimension, DoFHandler< dim > ::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< DoFHandler< dim > ::space_dimension > &data_postprocessor) |
void | add_data_vector (const DoFHandler< dim > &dof_handler, const VectorType &data, const DataPostprocessor< DoFHandler< dim > ::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 |
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 |
std::string | write_vtu_with_pvtu_record (const std::string &directory, const std::string &filename_without_extension, const unsigned int counter, const MPI_Comm &mpi_communicator, const unsigned int n_digits_for_counter=numbers::invalid_unsigned_int, const unsigned int n_groups=0) const |
void | write_svg (std::ostream &out) const |
void | write_deal_II_intermediate (std::ostream &out) const |
XDMFEntry | create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter, const std::string &h5_filename, const double cur_time, MPI_Comm comm) const |
XDMFEntry | create_xdmf_entry (const DataOutBase::DataOutFilter &data_filter, const std::string &h5_mesh_filename, const std::string &h5_solution_filename, const double cur_time, MPI_Comm comm) const |
void | write_xdmf_file (const std::vector< XDMFEntry > &entries, const std::string &filename, MPI_Comm comm) const |
void | write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter, const std::string &filename, MPI_Comm comm) const |
void | write_hdf5_parallel (const DataOutBase::DataOutFilter &data_filter, const bool write_mesh_file, const std::string &mesh_filename, const std::string &solution_filename, MPI_Comm comm) const |
void | write_filtered_data (DataOutBase::DataOutFilter &filtered_data) const |
void | write (std::ostream &out, const DataOutBase::OutputFormat output_format=DataOutBase::default_format) const |
void | set_default_format (const DataOutBase::OutputFormat default_format) |
void | set_flags (const FlagType &flags) |
std::string | default_suffix (const DataOutBase::OutputFormat output_format=DataOutBase::default_format) const |
void | parse_parameters (ParameterHandler &prm) |
std::size_t | memory_consumption () const |
Static Public Attributes | |
static constexpr unsigned int | spacedim = DoFHandlerType::space_dimension |
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) |
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 |
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.
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.
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 148 of file data_out.h.
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 166 of file data_out.h.
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 173 of file data_out.h.
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 181 of file data_out.h.
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 196 of file data_out.h.
Constructor.
Definition at line 69 of file data_out.cc.
|
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.
n_subdivisions | A 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. |
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 1071 of file data_out.cc.
|
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.
mapping
argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler. Definition at line 1083 of file data_out.cc.
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).
[in] | first_cell | A 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_cell | A 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.
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 1274 of file data_out.cc.
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:
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.
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 1288 of file data_out.cc.
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 1330 of file data_out.cc.
|
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 1339 of file data_out.cc.
|
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.
Definition at line 1348 of file data_out.cc.
|
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 1364 of file data_out.cc.
|
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 1381 of file data_out.cc.
|
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.
|
staticconstexpr |
Definition at line 157 of file data_out.h.
|
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 462 of file data_out.h.
|
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 471 of file data_out.h.