Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Types | Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
DataOutFaces< dim, DoFHandlerType > Class Template Reference

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

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

Public Types

using cell_iterator = typename DataOut_DoFData< DoFHandlerType, dimension - 1, dimension >::cell_iterator
 
using FaceDescriptor = typename std::pair< cell_iterator, unsigned int >
 
- Public Types inherited from DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension - 1, DoFHandler< dim > ::dimension >
enum  DataVectorType
 
using cell_iterator = typename Triangulation< DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension >::cell_iterator
 

Public Member Functions

 DataOutFaces (const bool surface_only=true)
 
virtual void build_patches (const unsigned int n_subdivisions=0)
 
virtual void build_patches (const Mapping< dimension > &mapping, const unsigned int n_subdivisions=0)
 
virtual FaceDescriptor first_face ()
 
virtual FaceDescriptor next_face (const FaceDescriptor &face)
 
- Public Member Functions inherited from DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension - 1, DoFHandler< dim > ::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 const unsigned int dimension = DoFHandlerType::dimension
 
static const unsigned int space_dimension = DoFHandlerType::space_dimension
 

Private Member Functions

void build_one_patch (const FaceDescriptor *cell_and_face, internal::DataOutFacesImplementation::ParallelData< dimension, dimension > &data, DataOutBase::Patch< dimension - 1, space_dimension > &patch)
 

Private Attributes

const bool surface_only
 

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< DoFHandler< dim >, DoFHandler< dim > ::dimension - 1, DoFHandler< dim > ::dimension >
using Patch = ::DataOutBase::Patch< patch_dim, patch_space_dim >
 
- Protected Member Functions inherited from DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension - 1, DoFHandler< dim > ::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< DoFHandler< dim > ::dimension, DoFHandler< dim > ::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 >
void validate_dataset_names () const
 
- Protected Attributes inherited from DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension - 1, DoFHandler< dim > ::dimension >
SmartPointer< const Triangulation< DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension > > triangulation
 
SmartPointer< const DoFHandler< dim > > dofs
 
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandler< dim > > > > dof_data
 
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandler< dim > > > > 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 DataOutFaces< dim, DoFHandlerType >

This class generates output from faces of a triangulation. It might be used to generate output only for the surface of the triangulation (this is the default of this class), or for all faces of active cells, as specified in the constructor. The output of this class is a set of patches (as defined by the class DataOutBase::Patch()), one for each face for which output is to be generated. These patches can then be written in several graphical data formats by the functions of the underlying classes.

Interface

The interface of this class is copied from the DataOut class. Furthermore, they share the common parent class DataOut_DoFData. See the reference of these two classes for a discussion of the interface.

Extending this class

The sequence of faces to generate patches from is generated in the same way as in the DataOut class; see there for a description of the respective interface. The functions generating the sequence of faces which shall be used to generate output, are called first_face() and next_face().

Since we need to initialize objects of type FEValues with the faces generated from these functions, it is not sufficient that they only return face iterators. Rather, we need a pair of cell and the number of the face, as the values of finite element fields need not necessarily be unique on a face (think of discontinuous finite elements, where the value of the finite element field depend on the direction from which you approach a face, thus it is necessary to use a pair of cell and face, rather than only a face iterator). Therefore, this class defines an alias which creates a type FaceDescriptor that is an abbreviation for a pair of cell iterator and face number. The functions first_face and next_face operate on objects of this type.

Extending this class might, for example, be useful if you only want output from certain portions of the boundary, e.g. as indicated by the boundary indicator of the respective faces. However, it is also conceivable that one generates patches not from boundary faces, but from interior faces that are selected due to other criteria; one application might be to use only those faces where one component of the solution attains a certain value, in order to display the values of other solution components on these faces. Other applications certainly exist, for which the author is not imaginative enough.

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.
Todo:
Reimplement this whole class using actual FEFaceValues and MeshWorker.
Author
Wolfgang Bangerth, Guido Kanschat, 2000, 2011

Definition at line 117 of file data_out_faces.h.

Member Typedef Documentation

◆ cell_iterator

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
using DataOutFaces< dim, DoFHandlerType >::cell_iterator = typename DataOut_DoFData<DoFHandlerType, dimension - 1, dimension>::cell_iterator

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

Definition at line 140 of file data_out_faces.h.

◆ FaceDescriptor

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
using DataOutFaces< dim, DoFHandlerType >::FaceDescriptor = typename std::pair<cell_iterator, unsigned int>

Declare a way to describe a face which we would like to generate output for. The usual way would, of course, be to use an object of type DoFHandler<dim>::face_iterator, but since we have to describe faces to objects of type FEValues, we can only represent faces by pairs of a cell and the number of the face. This pair is here aliased to a name that is better to type.

Definition at line 197 of file data_out_faces.h.

Constructor & Destructor Documentation

◆ DataOutFaces()

template<int dim, typename DoFHandlerType >
DataOutFaces< dim, DoFHandlerType >::DataOutFaces ( const bool  surface_only = true)

Constructor determining whether a surface mesh (default) or the whole wire basket is written.

Definition at line 84 of file data_out_faces.cc.

Member Function Documentation

◆ build_patches() [1/2]

template<int dim, typename DoFHandlerType >
void DataOutFaces< 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 face 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_subdivisionsSee DataOut::build_patches() for an extensive description of this parameter.

Definition at line 321 of file data_out_faces.cc.

◆ build_patches() [2/2]

template<int dim, typename DoFHandlerType >
void DataOutFaces< dim, DoFHandlerType >::build_patches ( const Mapping< dimension > &  mapping,
const unsigned int  n_subdivisions = 0 
)
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.

Even for non-curved cells the mapping argument can be used for the 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 331 of file data_out_faces.cc.

◆ first_face()

template<int dim, typename DoFHandlerType >
DataOutFaces< dim, DoFHandlerType >::FaceDescriptor DataOutFaces< dim, DoFHandlerType >::first_face
virtual

Return the first face which we want output for. The default implementation returns the first face of a (locally owned) active cell or, if the surface_only option was set in the destructor (as is the default), the first such face that is located on the boundary.

If you want to use a different logic to determine which faces should contribute to the creation of graphical output, you can overload this function in a derived class.

Definition at line 420 of file data_out_faces.cc.

◆ next_face()

template<int dim, typename DoFHandlerType >
DataOutFaces< dim, DoFHandlerType >::FaceDescriptor DataOutFaces< dim, DoFHandlerType >::next_face ( const FaceDescriptor face)
virtual

Return the next face after which we want output for. If there are no more faces, dofs->end() is returned as the first component of the return value.

The default implementation returns the next face of a (locally owned) active cell, or the next such on the boundary (depending on whether the surface_only option was provided to the constructor).

This function traverses the mesh active cell by active cell (restricted to locally owned cells), and then through all faces of the cell. As a result, interior faces are output twice, a feature that is useful for discontinuous Galerkin methods or if a DataPostprocessor is used that might produce results that are discontinuous between cells).

This function can be overloaded in a derived class to select a different set of faces. Note that the default implementation assumes that the given face is active, which is guaranteed as long as first_face() is also used from the default implementation. Overloading only one of the two functions should be done with care.

Definition at line 441 of file data_out_faces.cc.

◆ build_one_patch()

template<int dim, typename DoFHandlerType >
void DataOutFaces< dim, DoFHandlerType >::build_one_patch ( const FaceDescriptor cell_and_face,
internal::DataOutFacesImplementation::ParallelData< dimension, dimension > &  data,
DataOutBase::Patch< dimension - 1, space_dimension > &  patch 
)
private

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

Definition at line 94 of file data_out_faces.cc.

Member Data Documentation

◆ dimension

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
const unsigned int DataOutFaces< dim, DoFHandlerType >::dimension = DoFHandlerType::dimension
static

An abbreviation for the dimension of the DoFHandler object we work with. Faces are then dimension-1 dimensional objects.

Definition at line 126 of file data_out_faces.h.

◆ space_dimension

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
const unsigned int DataOutFaces< dim, DoFHandlerType >::space_dimension = DoFHandlerType::space_dimension
static

An abbreviation for the spatial dimension within which the triangulation and DoFHandler are embedded in.

Definition at line 132 of file data_out_faces.h.

◆ surface_only

template<int dim, typename DoFHandlerType = DoFHandler<dim>>
const bool DataOutFaces< dim, DoFHandlerType >::surface_only
private

Parameter deciding between surface meshes and full wire basket.

Definition at line 241 of file data_out_faces.h.


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