
Public Types | |
| typedef DataOut_DoFData< DH, DH::dimension-1, DH::dimension > ::cell_iterator | cell_iterator |
| typedef std::pair < cell_iterator, unsigned int > | FaceDescriptor |
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) |
Static Public Member Functions | |
| ::ExceptionBase & | ExcInvalidNumberOfSubdivisions (int arg1) throw (errortext << "The number of subdivisions per patch, " << arg1 << ", is not valid." ) |
| ::ExceptionBase & | ExcCellNotActiveForCellData () |
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 |
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.
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.
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.
Definition at line 119 of file data_out_faces.h.
| 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.
Reimplemented from DataOut_DoFData< DH, DH::dimension-1, DH::dimension >.
Definition at line 129 of file data_out_faces.h.
| 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 220 of file data_out_faces.h.
| DataOutFaces< dim, DH >::DataOutFaces | ( | const bool | surface_only = true ) |
Constructor determining whether a surface mesh (default) or the whole wire basket is written.
| 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. The default number of threads to be used to build the patches is set to multithread_info.n_default_threads.
| 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.
mapping argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler. | 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.
| 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.
| ::ExceptionBase& DataOutFaces< dim, DH >::ExcInvalidNumberOfSubdivisions | ( | int | arg1 ) | throw (errortext << "The number of subdivisions per patch, " << arg1 << ", is not valid." ) [static] |
Exception
| ::ExceptionBase& DataOutFaces< dim, DH >::ExcCellNotActiveForCellData | ( | ) | [static] |
Exception
| 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.
const bool DataOutFaces< dim, DH >::surface_only [private] |
Parameter deciding between surface meshes and full wire basket.
Definition at line 291 of file data_out_faces.h.
documentation generated on Fri Feb 3 2012 06:04:06 by
doxygen
1.7.2