Loading [MathJax]/extensions/TeX/newcommand.js
 Reference documentation for deal.II version 9.4.1
\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\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
Particles::DataOut< dim, spacedim > Class Template Reference

#include <deal.II/particles/data_out.h>

Inheritance diagram for Particles::DataOut< dim, spacedim >:
[legend]

Public Member Functions

 DataOut ()=default
 
 ~DataOut ()=default
 
void build_patches (const Particles::ParticleHandler< dim, spacedim > &particles, const std::vector< std::string > &data_component_names={}, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretations={})
 
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)
 
std::size_t memory_consumption () const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Protected Member Functions

virtual const std::vector< DataOutBase::Patch< 0, spacedim > > & get_patches () const override
 
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
 
void validate_dataset_names () const
 

Protected Attributes

unsigned int default_subdivisions
 

Private Attributes

std::vector< DataOutBase::Patch< 0, spacedim > > patches
 
std::vector< std::string > dataset_names
 
std::vector< DataComponentInterpretation::DataComponentInterpretationdata_component_interpretations
 
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
 

Detailed Description

template<int dim, int spacedim = dim>
class Particles::DataOut< dim, spacedim >

This class generates graphical output for the particles stored by a ParticleHandler object. From a particle handler, it generates patches which can then be used to write traditional output files. This class currently only supports witing the particle position and their ID and does not allow to write the properties attached to the particles

Definition at line 44 of file data_out.h.

Constructor & Destructor Documentation

◆ DataOut()

template<int dim, int spacedim = dim>
Particles::DataOut< dim, spacedim >::DataOut ( )
default

Default constructor for the Particles::DataOut class.

◆ ~DataOut()

template<int dim, int spacedim = dim>
Particles::DataOut< dim, spacedim >::~DataOut ( )
default

Destructor for the Particles::DataOut class.

Member Function Documentation

◆ build_patches()

template<int dim, int spacedim>
void DataOut< dim, spacedim >::build_patches ( const Particles::ParticleHandler< dim, spacedim > &  particles,
const std::vector< std::string > &  data_component_names = {},
const std::vector< DataComponentInterpretation::DataComponentInterpretation > &  data_component_interpretations = {} 
)

Build the patches for a given particle handler.

Parameters
[in]particlesA particle handler for which the patches will be built. A dim=0 patch is built for each particle. The position of the particle is used to build the node position and the ID of the particle is added as a single data element.
[in]data_component_namesAn optional vector of strings that describe the properties of each particle. Particle properties will only be written if this vector is provided.
[in]data_component_interpretationsAn optional vector that controls if the particle properties are interpreted as scalars, vectors, or tensors. Has to be of the same length as data_component_names.

Definition at line 28 of file data_out.cc.

◆ get_patches()

template<int dim, int spacedim>
const std::vector< DataOutBase::Patch< 0, spacedim > > & DataOut< dim, spacedim >::get_patches
overrideprotectedvirtual

Returns the patches built by the data_out class which was previously built using a particle handler

Implements DataOutInterface< dim, spacedim >.

Definition at line 107 of file data_out.cc.

◆ get_dataset_names()

template<int dim, int spacedim>
std::vector< std::string > DataOut< dim, spacedim >::get_dataset_names
overrideprotectedvirtual

Virtual function through which the names of data sets are obtained from this class

Implements DataOutInterface< dim, spacedim >.

Definition at line 116 of file data_out.cc.

◆ get_nonscalar_data_ranges()

template<int dim, int spacedim>
std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > DataOut< dim, spacedim >::get_nonscalar_data_ranges
overrideprotectedvirtual

Overload of the respective DataOutInterface::get_nonscalar_data_ranges() function. See there for a more extensive documentation. This function is a reimplementation of the function DataOut_DoFData::get_nonscalar_data_ranges().

Reimplemented from DataOutInterface< dim, spacedim >.

Definition at line 129 of file data_out.cc.

Member Data Documentation

◆ patches

template<int dim, int spacedim = dim>
std::vector<DataOutBase::Patch<0, spacedim> > Particles::DataOut< dim, spacedim >::patches
private

This is a vector 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 115 of file data_out.h.

◆ dataset_names

template<int dim, int spacedim = dim>
std::vector<std::string> Particles::DataOut< dim, spacedim >::dataset_names
private

A vector of field names for all data components stored in patches.

Definition at line 120 of file data_out.h.

◆ data_component_interpretations

template<int dim, int spacedim = dim>
std::vector<DataComponentInterpretation::DataComponentInterpretation> Particles::DataOut< dim, spacedim >::data_component_interpretations
private

A vector that for each of the data components of the current data set indicates whether they are scalar fields, parts of a vector-field, or any of the other supported kinds of data.

Definition at line 128 of file data_out.h.


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