Reference documentation for deal.II version 9.0.0
|
#include <deal.II/numerics/data_out_rotation.h>
Public Types | |
typedef DataOut_DoFData< DoFHandlerType, dimension+1 >::cell_iterator | cell_iterator |
Public Types inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension+1 > | |
enum | DataVectorType |
typedef Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension >::cell_iterator | cell_iterator |
Public Member Functions | |
virtual void | build_patches (const unsigned int n_patches_per_circle, const unsigned int n_subdivisions=0) |
virtual cell_iterator | first_cell () |
virtual cell_iterator | next_cell (const cell_iterator &cell) |
Public Member Functions inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension+1 > | |
DataOut_DoFData () | |
virtual | ~DataOut_DoFData () |
void | attach_dof_handler (const DoFHandlerType &) |
void | attach_triangulation (const Triangulation< DoFHandlerType::dimension, DoFHandlerType::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 DoFHandlerType &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 DoFHandlerType &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< DoFHandlerType::space_dimension > &data_postprocessor) |
void | add_data_vector (const DoFHandlerType &dof_handler, const VectorType &data, const DataPostprocessor< DoFHandlerType::space_dimension > &data_postprocessor) |
void | clear_data_vectors () |
void | clear_input_data_references () |
void | merge_patches (const DataOut_DoFData< DoFHandlerType2, patch_dim, patch_dim > &source, const Point< patch_dim > &shift=Point< patch_dim >()) |
virtual void | clear () |
std::size_t | memory_consumption () const |
Public Member Functions inherited from DataOutInterface< patch_dim, patch_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 char *filename, MPI_Comm comm) const |
void | write_pvtu_record (std::ostream &out, const 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 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 Member Functions | |
static ::ExceptionBase & | ExcRadialVariableHasNegativeValues (double arg1) |
Static Public Member Functions inherited from DataOutInterface< patch_dim, patch_dim > | |
static void | declare_parameters (ParameterHandler &prm) |
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 cell_iterator *cell, internal::DataOutRotationImplementation::ParallelData< dimension, space_dimension > &data, std::vector< DataOutBase::Patch< dimension+1, space_dimension+1 > > &my_patches) |
Additional Inherited Members | |
Protected Types inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension+1 > | |
typedef ::DataOutBase::Patch< patch_dim, patch_dim > | Patch |
Protected Member Functions inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension+1 > | |
virtual const std::vector< Patch > & | get_patches () const |
virtual std::vector< std::string > | get_dataset_names () const |
std::vector< std::shared_ptr<::hp::FECollection< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > > | get_fes () const |
virtual std::vector< std::tuple< unsigned int, unsigned int, std::string > > | get_vector_data_ranges () const |
Protected Member Functions inherited from DataOutInterface< patch_dim, patch_dim > | |
void | validate_dataset_names () const |
Protected Attributes inherited from DataOut_DoFData< DoFHandlerType, DoFHandlerType::dimension+1 > | |
SmartPointer< const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > | triangulation |
SmartPointer< const DoFHandlerType > | dofs |
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType > > > | dof_data |
std::vector< std::shared_ptr< internal::DataOutImplementation::DataEntryBase< DoFHandlerType > > > | cell_data |
std::vector< Patch > | patches |
Protected Attributes inherited from DataOutInterface< patch_dim, patch_dim > | |
unsigned int | default_subdivisions |
This class generates output in the full domain of computations that were done using rotational symmetry of domain and solution. In particular, if a computation of a three dimensional problem with rotational symmetry around the z-axis
(i.e. in the r-z-plane
) was done, then this class can be used to generate the output in the original x-y-z
space. In order to do so, it generates from each cell in the computational mesh a cell in the space with dimension one greater than that of the DoFHandler object. The resulting output will then consist of hexahedra forming an object that has rotational symmetry around the z-axis. As most graphical programs can not represent ring-like structures, the angular (rotation) variable is discretized into a finite number of intervals as well; the number of these intervals must be given to the build_patches
function. It is noted, however, that while this function generates nice pictures of the whole domain, it often produces very large output files.
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 and how to extend it by deriving further classes from this class.
The one coordinate in the triangulation used by the DoFHandler object passed to this class is taken as the radial variable, and the output will then be either a circle or a ring domain. It is in the user's responsibility to assure that the radial coordinate only attains non- negative values.
We consider the computation (represented by the DoFHandler object that is attached to this class) to have happened in the r-z-plane
, where r
is the radial variable and z
denotes the axis of revolution around which the solution is symmetric. The output is in x-y-z
space, where the radial dependence is transformed to the x-y
plane. At present, it is not possible to exchange the meaning of the first and second variable of the plane in which the simulation was made, i.e. generate output from a simulation where the first variable denoted the symmetry axis, and the second denoted the radial variable. You have to take that into account when first programming your application.
It is in the responsibility of the user to make sure that the radial variable attains only non-negative values.
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.Definition at line 119 of file data_out_rotation.h.
typedef DataOut_DoFData<DoFHandlerType,dimension+1>::cell_iterator DataOutRotation< dim, DoFHandlerType >::cell_iterator |
Typedef to the iterator type of the dof handler class under consideration.
Definition at line 138 of file data_out_rotation.h.
|
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_patches_per_circle | Denotes into how many intervals the angular (rotation) variable is to be subdivided. |
n_subdivisions | See DataOut::build_patches() for an extensive description of this parameter. |
Definition at line 412 of file data_out_rotation.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 500 of file data_out_rotation.cc.
|
virtual |
Return the next cell after cell
which we want output for. If there are no more cells, dofs->end()
shall be returned.
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 508 of file data_out_rotation.cc.
|
private |
Build all of the patches that correspond to the cell given in the first argument. Use the second argument as scratch space for parallel invocation in WorkStream, and put the results into the last argument.
Definition at line 96 of file data_out_rotation.cc.
|
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_rotation.h.
|
static |
An abbreviation for the spatial dimension within which the triangulation and DoFHandler are embedded in.
Definition at line 132 of file data_out_rotation.h.