Reference documentation for deal.II version 8.1.0
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
DataOutFaces< dim, DH > Class Template Reference

#include <data_out_faces.h>

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

Public Types

typedef DataOut_DoFData< DH, DH::dimension-1, DH::dimension >::cell_iterator cell_iterator
 
typedef std::pair< cell_iterator, unsigned intFaceDescriptor
 
- Public Types inherited from DataOut_DoFData< DH, DH::dimension-1, DH::dimension >
enum  DataVectorType
 
typedef Triangulation< DH::dimension, DH::space_dimension >::cell_iterator cell_iterator
 
typedef Triangulation< DH::dimension, DH::space_dimension >::active_cell_iterator active_cell_iterator
 
- Public Types inherited from DataOutInterface< patch_dim, patch_space_dim >
enum  OutputFormat
 

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< DH::dimension > &mapping, const unsigned int n_subdivisions=0)
 
virtual FaceDescriptor first_face ()
 
virtual FaceDescriptor next_face (const FaceDescriptor &face)
 
 DeclException1 (ExcInvalidNumberOfSubdivisions, int,<< "The number of subdivisions per patch, "<< arg1<< ", is not valid.")
 
 DeclException0 (ExcCellNotActiveForCellData)
 
- Public Member Functions inherited from DataOut_DoFData< DH, DH::dimension-1, DH::dimension >
 DataOut_DoFData ()
 
virtual ~DataOut_DoFData ()
 
void attach_dof_handler (const DH &)
 
void attach_triangulation (const Triangulation< DH::dimension, DH::space_dimension > &)
 
void add_data_vector (const VECTOR &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 VECTOR &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 DH &dof_handler, const VECTOR &data, const std::vector< std::string > &names, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_data_vector (const DH &dof_handler, const VECTOR &data, const std::string &name, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
 
void add_data_vector (const VECTOR &data, const DataPostprocessor< DH::space_dimension > &data_postprocessor)
 
void add_data_vector (const DH &dof_handler, const VECTOR &data, const DataPostprocessor< DH::space_dimension > &data_postprocessor)
 
void clear_data_vectors ()
 
void clear_input_data_references ()
 
void merge_patches (const DataOut_DoFData< DH2, 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
 
 DeclException0 (ExcNoTriangulationSelected)
 
 DeclException0 (ExcNoDoFHandlerSelected)
 
 DeclException0 (ExcDataPostprocessingIsNotPossibleForCellData)
 
 DeclException0 (ExcOldDataStillPresent)
 
 DeclException0 (ExcNoPatches)
 
 DeclException0 (ExcIncompatibleDatasetNames)
 
 DeclException0 (ExcIncompatiblePatchLists)
 
 DeclException3 (ExcInvalidVectorSize, int, int, int,<< "The vector has size "<< arg1<< " but the DoFHandler objects says there are "<< arg2<< " degrees of freedom and there are "<< arg3<< " active cells.")
 
 DeclException2 (ExcInvalidCharacter, std::string, size_t,<< "Please use only the characters [a-zA-Z0-9_<>()] for"<< std::endl<< "description strings since some graphics formats will only accept these."<< std::endl<< "The string you gave was <"<< arg1<< ">, the invalid character is <"<< arg1[arg2]<< ">."<< std::endl)
 
 DeclException2 (ExcInvalidNumberOfNames, int, int,<< "You have to give one name per component in your "<< "data vector. The number you gave was "<< arg1<< ", but the number of components is "<< arg2)
 
 DeclException2 (ExcInvalidVectorDeclaration, int, std::string,<< "When declaring that a number of components in a data\n"<< "set to be output logically form a vector instead of\n"<< "simply a set of scalar fields, you need to specify\n"<< "this for all relevant components. Furthermore,\n"<< "vectors must always consist of exactly <dim>\n"<< "components. However, the vector component at\n"<< "position "<< arg1<< " with name <"<< arg2<< "> does not satisfy these conditions.")
 
- Public Member Functions inherited from DataOutInterface< patch_dim, patch_space_dim >
 DataOutInterface ()
 
virtual ~DataOutInterface ()
 
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 char *filename, MPI_Comm comm) const
 
void write_pvtu_record (std::ostream &out, const std::vector< std::string > &piece_names) const
 
void write_pvd_record (std::ostream &out, const std::vector< std::pair< double, std::string > > &times_and_names) const
 
void write_visit_record (std::ostream &out, const std::vector< std::string > &piece_names) const
 
void write_visit_record (std::ostream &out, const std::vector< 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 std::string &h5_filename, const double cur_time, MPI_Comm comm) const DEAL_II_DEPRECATED
 
XDMFEntry create_xdmf_entry (const DataOutFilter &data_filter, const std::string &h5_filename, const double cur_time, MPI_Comm comm) const
 
XDMFEntry create_xdmf_entry (const 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 std::string &filename, MPI_Comm comm) const DEAL_II_DEPRECATED
 
void write_hdf5_parallel (const DataOutFilter &data_filter, const std::string &filename, MPI_Comm comm) const
 
void write_hdf5_parallel (const 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 (DataOutFilter &filtered_data) const
 
void write (std::ostream &out, const OutputFormat output_format=default_format) const
 
void set_default_format (const OutputFormat default_format)
 
void set_flags (const DXFlags &dx_flags)
 
void set_flags (const UcdFlags &ucd_flags)
 
void set_flags (const GnuplotFlags &gnuplot_flags)
 
void set_flags (const PovrayFlags &povray_flags)
 
void set_flags (const EpsFlags &eps_flags)
 
void set_flags (const GmvFlags &gmv_flags)
 
void set_flags (const TecplotFlags &tecplot_flags)
 
void set_flags (const VtkFlags &vtk_flags)
 
void set_flags (const SvgFlags &svg_flags)
 
void set_flags (const Deal_II_IntermediateFlags &deal_II_intermediate_flags)
 
std::string default_suffix (const OutputFormat output_format=default_format) const
 
void parse_parameters (ParameterHandler &prm)
 
std::size_t memory_consumption () const
 

Private Member Functions

void build_one_patch (const FaceDescriptor *cell_and_face, internal::DataOutFaces::ParallelData< DH::dimension, DH::dimension > &data, DataOutBase::Patch< DH::dimension-1, DH::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< DH, DH::dimension-1, DH::dimension >
typedef ::DataOutBase::Patch< patch_dim, patch_space_dim > Patch
 
- Protected Member Functions inherited from DataOut_DoFData< DH, DH::dimension-1, DH::dimension >
virtual const std::vector< Patch > & get_patches () const
 
virtual std::vector< std::string > get_dataset_names () const
 
std::vector< std_cxx1x::shared_ptr<::hp::FECollection< DH::dimension, DH::space_dimension > > > get_finite_elements () const
 
virtual std::vector< std_cxx1x::tuple< unsigned int, unsigned int, std::string > > get_vector_data_ranges () const
 
- Protected Attributes inherited from DataOut_DoFData< DH, DH::dimension-1, DH::dimension >
SmartPointer< const Triangulation< DH::dimension, DH::space_dimension > > triangulation
 
SmartPointer< const DH > dofs
 
std::vector< std_cxx1x::shared_ptr< internal::DataOut::DataEntryBase< DH > > > dof_data
 
std::vector< std_cxx1x::shared_ptr< internal::DataOut::DataEntryBase< DH > > > cell_data
 
std::vector< Patchpatches
 
- Protected Attributes inherited from DataOutInterface< patch_dim, patch_space_dim >
unsigned int default_subdivisions
 

Detailed Description

template<int dim, class DH = DoFHandler<dim>>
class DataOutFaces< dim, DH >

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 a typedef 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.

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

template<int dim, class DH = DoFHandler<dim>>
typedef DataOut_DoFData<DH,DH::dimension-1, DH::dimension>::cell_iterator DataOutFaces< dim, DH >::cell_iterator

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

Definition at line 127 of file data_out_faces.h.

template<int dim, class DH = DoFHandler<dim>>
typedef std::pair<cell_iterator,unsigned int> DataOutFaces< dim, DH >::FaceDescriptor

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 215 of file data_out_faces.h.

Constructor & Destructor Documentation

template<int dim, class DH = DoFHandler<dim>>
DataOutFaces< dim, DH >::DataOutFaces ( const bool  surface_only = true)

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

Member Function Documentation

template<int dim, class DH = DoFHandler<dim>>
virtual void DataOutFaces< dim, DH >::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. See the general documentation of this class for further information.

The function supports multithreading, if deal.II is compiled in multithreading mode.

template<int dim, class DH = DoFHandler<dim>>
virtual void DataOutFaces< dim, DH >::build_patches ( const Mapping< DH::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.
template<int dim, class DH = DoFHandler<dim>>
virtual FaceDescriptor DataOutFaces< dim, DH >::first_face ( )
virtual

Return the first face which we want output for. The default implementation returns the first face of an active cell or the first such on the boundary.

For more general sets, overload this function in a derived class.

template<int dim, class DH = DoFHandler<dim>>
virtual FaceDescriptor DataOutFaces< dim, DH >::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 an active cell, or the next such on the boundary.

This function traverses the mesh cell by cell (active only), and then through all faces of the cell. As a result, interior faces are output twice.

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.

template<int dim, class DH = DoFHandler<dim>>
DataOutFaces< dim, DH >::DeclException1 ( ExcInvalidNumberOfSubdivisions  ,
int  ,
<< "The number of subdivisions per  patch,
"<< arg1<< "  ,
is not valid."   
)

Exception

template<int dim, class DH = DoFHandler<dim>>
DataOutFaces< dim, DH >::DeclException0 ( ExcCellNotActiveForCellData  )

Exception

template<int dim, class DH = DoFHandler<dim>>
void DataOutFaces< dim, DH >::build_one_patch ( const FaceDescriptor cell_and_face,
internal::DataOutFaces< dim, DH >::ParallelData< DH::dimension, DH::dimension > &  data,
DataOutBase::Patch< DH::dimension-1, DH::space_dimension > &  patch 
)
private

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

Member Data Documentation

template<int dim, class DH = DoFHandler<dim>>
const bool DataOutFaces< dim, DH >::surface_only
private

Parameter deciding between surface meshes and full wire basket.

Definition at line 286 of file data_out_faces.h.


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