Reference documentation for deal.II version 9.4.1
|
#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< dim, dim, spacedim, spacedim >::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 { type_dof_data , type_cell_data , type_automatic } |
Public Member Functions | |
DataOut () | |
virtual void | build_patches (const unsigned int n_subdivisions=0) |
virtual void | build_patches (const Mapping< dim, spacedim > &mapping, const unsigned int n_subdivisions=0, const CurvedCellRegion curved_region=curved_boundary) |
virtual void | build_patches (const hp::MappingCollection< dim, spacedim > &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 |
void | attach_dof_handler (const DoFHandler< dim, spacedim > &) |
void | attach_triangulation (const Triangulation< dim, spacedim > &) |
template<class VectorType > | |
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={}) |
template<class VectorType > | |
void | add_data_vector (const VectorType &data, const std::string &name, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={}) |
template<class VectorType > | |
void | add_data_vector (const DoFHandler< dim, spacedim > &dof_handler, const VectorType &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={}) |
template<class VectorType > | |
void | add_data_vector (const DoFHandler< dim, spacedim > &dof_handler, const VectorType &data, const std::string &name, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={}) |
template<class VectorType > | |
void | add_data_vector (const VectorType &data, const DataPostprocessor< spacedim > &data_postprocessor) |
template<class VectorType > | |
void | add_data_vector (const DoFHandler< dim, spacedim > &dof_handler, const VectorType &data, const DataPostprocessor< spacedim > &data_postprocessor) |
template<class VectorType > | |
void | add_mg_data_vector (const DoFHandler< dim, spacedim > &dof_handler, const MGLevelObject< VectorType > &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >()) |
template<class VectorType > | |
void | add_mg_data_vector (const DoFHandler< dim, spacedim > &dof_handler, const MGLevelObject< VectorType > &data, const std::string &name) |
void | clear_data_vectors () |
void | clear_input_data_references () |
template<int dim2, int spacedim2> | |
void | merge_patches (const DataOut_DoFData< dim2, patch_dim, spacedim2, patch_spacedim > &source, const Point< patch_spacedim > &shift=Point< patch_spacedim >()) |
template<typename DoFHandlerType2 > | |
void | merge_patches (const DataOut_DoFData< DoFHandlerType2::dimension, patch_dim, DoFHandlerType2::space_dimension, patch_spacedim > &source, const Point< patch_spacedim > &shift=Point< patch_spacedim >()) |
virtual void | clear () |
std::size_t | memory_consumption () const |
virtual const std::vector< Patch > & | get_patches () const override |
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) |
template<typename FlagType > | |
void | set_flags (const FlagType &flags) |
std::string | default_suffix (const DataOutBase::OutputFormat output_format=DataOutBase::default_format) const |
void | parse_parameters (ParameterHandler &prm) |
Static Public Member Functions | |
static void | declare_parameters (ParameterHandler &prm) |
Protected Types | |
using | Patch = ::DataOutBase::Patch< patch_dim, patch_spacedim > |
Protected Member Functions | |
virtual std::vector< std::string > | get_dataset_names () const override |
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > | get_nonscalar_data_ranges () const override |
std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > | get_fes () const |
void | validate_dataset_names () const |
Protected Attributes | |
SmartPointer< const Triangulation< dim, spacedim > > | triangulation |
SmartPointer< const DoFHandler< dim, spacedim > > | dofs |
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< dim, spacedim > > > | dof_data |
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< dim, spacedim > > > | cell_data |
std::vector< Patch > | patches |
unsigned int | default_subdivisions |
Private Member Functions | |
void | build_one_patch (const std::pair< cell_iterator, unsigned int > *cell_and_index, internal::DataOutImplementation::ParallelData< dim, spacedim > &scratch_data, const unsigned int n_subdivisions, const CurvedCellRegion curved_cell_region) |
template<class VectorType > | |
void | add_data_vector_internal (const DoFHandler< dim, spacedim > *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 |
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 private std::function objects first_cell_function() and next_cell_function(). 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.
using DataOut< dim, spacedim >::cell_iterator = typename DataOut_DoFData<dim, dim, spacedim, spacedim>::cell_iterator |
Typedef to the iterator type of the dof handler class under consideration.
Definition at line 154 of file data_out.h.
using DataOut< dim, spacedim >::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 161 of file data_out.h.
using DataOut< dim, spacedim >::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 168 of file data_out.h.
|
protectedinherited |
Abbreviate the somewhat lengthy name for the Patch class.
Definition at line 963 of file data_out_dof_data.h.
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 185 of file data_out.h.
|
inherited |
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).
Enumerator | |
---|---|
type_dof_data | Data vector entries are associated to degrees of freedom |
type_cell_data | Data vector entries are one per grid cell |
type_automatic | Find out automatically |
Definition at line 615 of file data_out_dof_data.h.
Constructor.
Definition at line 69 of file data_out.cc.
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 1064 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.
Definition at line 1078 of file data_out.cc.
|
virtual |
Same as above, but for hp::MappingCollection.
Definition at line 1091 of file data_out.cc.
void DataOut< dim, spacedim >::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 1275 of file data_out.cc.
void DataOut< dim, spacedim >::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 1289 of file data_out.cc.
const std::pair< typename DataOut< dim, spacedim >::FirstCellFunctionType, typename DataOut< dim, spacedim >::NextCellFunctionType > DataOut< dim, spacedim >::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 1331 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 97 of file data_out.cc.
|
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.
|
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.
|
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.
Definition at line 1082 of file data_out_dof_data.h.
|
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 1063 of file data_out_dof_data.h.
|
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 1121 of file data_out_dof_data.h.
|
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 1100 of file data_out_dof_data.h.
|
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.
Definition at line 1141 of file data_out_dof_data.h.
|
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.
|
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.
|
inherited |
Scalar version of the function above.
|
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.
|
inherited |
Release pointers to all input data elements, i.e. pointers to 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 must not be deleted before the output thread is finished.
|
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 1155 of file data_out_dof_data.h.
|
inherited |
Definition at line 1246 of file data_out_dof_data.h.
|
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.
|
inherited |
Determine an estimate for the memory consumption (in bytes) of this object.
|
overridevirtualinherited |
Function by which the base class's functions get to know what patches they shall write to a file.
Implements DataOutInterface< dim, spacedim >.
Reimplemented in DataOutResample< dim, patch_dim, spacedim >.
|
overrideprotectedvirtualinherited |
Virtual function through which the names of data sets are obtained by the output functions of the base class.
Implements DataOutInterface< dim, spacedim >.
|
overrideprotectedvirtualinherited |
Overload of the respective DataOutInterface::get_nonscalar_data_ranges() function. See there for a more extensive documentation.
Reimplemented from DataOutInterface< dim, spacedim >.
|
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.
|
privateinherited |
Common function called by the four public add_data_vector methods.
|
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 426 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 435 of file data_out.h.
|
protectedinherited |
Pointer to the triangulation object.
Definition at line 968 of file data_out_dof_data.h.
|
protectedinherited |
Pointer to the optional handler object.
Definition at line 973 of file data_out_dof_data.h.
|
protectedinherited |
List of data elements with vectors of values for each degree of freedom.
Definition at line 980 of file data_out_dof_data.h.
|
protectedinherited |
List of data elements with vectors of values for each cell.
Definition at line 987 of file data_out_dof_data.h.
|
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 994 of file data_out_dof_data.h.