dune-vtk  0.2
Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
Dune::VtkWriterInterface< GV, DC > Class Template Referenceabstract

Interface for file writers for the Vtk XML file formats. More...

#include <dune/vtk/vtkwriterinterface.hh>

Inheritance diagram for Dune::VtkWriterInterface< GV, DC >:
Inheritance graph

Public Types

using GridView = GV
 
using DataCollector = DC
 

Public Member Functions

 VtkWriterInterface (GridView const &gridView, Vtk::FormatTypes format=Vtk::FormatTypes::BINARY, Vtk::DataTypes datatype=Vtk::DataTypes::FLOAT32, Vtk::DataTypes headertype=Vtk::DataTypes::UINT32)
 Constructor, passes the gridView to the DataCollector. More...
 
 VtkWriterInterface (DataCollector &dataCollector, Vtk::FormatTypes format=Vtk::FormatTypes::BINARY, Vtk::DataTypes datatype=Vtk::DataTypes::FLOAT32, Vtk::DataTypes headertype=Vtk::DataTypes::UINT32)
 Constructor, wraps the passed DataCollector in a non-destroying shared_ptr. More...
 
 VtkWriterInterface (std::shared_ptr< DataCollector > dataCollector, Vtk::FormatTypes format=Vtk::FormatTypes::BINARY, Vtk::DataTypes datatype=Vtk::DataTypes::FLOAT32, Vtk::DataTypes headertype=Vtk::DataTypes::UINT32)
 Constructor, stores the passed DataCollector. More...
 
virtual std::string write (std::string const &fn, std::optional< std::string > dir={}) const override
 Write the attached data to the file. More...
 
template<class Function , class... Args>
VtkWriterInterfaceaddPointData (Function &&fct, Args &&... args)
 Attach point data to the writer. More...
 
template<class Function , class... Args>
VtkWriterInterfaceaddCellData (Function &&fct, Args &&... args)
 Attach cell data to the writer. More...
 
void setFormat (Vtk::FormatTypes format)
 
void setDatatype (Vtk::DataTypes datatype)
 Sets the global datatype used for coordinates and other global float values. More...
 
void setHeadertype (Vtk::DataTypes datatype)
 Sets the integer type used in binary data headers. More...
 
void setCompressor (Vtk::CompressorTypes compressor, int level=-1)
 

Protected Types

enum  PositionTypes { POINT_DATA , CELL_DATA }
 
using VtkFunction = Dune::Vtk::Function< GridView >
 
using pos_type = typename std::ostream::pos_type
 

Protected Member Functions

void writeData (std::ofstream &out, std::vector< pos_type > &offsets, VtkFunction const &fct, PositionTypes type, std::optional< std::size_t > timestep={}) const
 
void writeDataAppended (std::ofstream &out, std::vector< std::uint64_t > &blocks) const
 
void writePoints (std::ofstream &out, std::vector< pos_type > &offsets, std::optional< std::size_t > timestep={}) const
 
void writeAppended (std::ofstream &out, std::vector< pos_type > const &offsets) const
 
template<class HeaderType , class FloatType >
std::uint64_t writeValuesAppended (std::ofstream &out, std::vector< FloatType > const &values) const
 
template<class T >
void writeValuesAscii (std::ofstream &out, std::vector< T > const &values) const
 
void writeHeader (std::ofstream &out, std::string const &type) const
 
std::string getNames (std::vector< VtkFunction > const &data) const
 Return PointData/CellData attributes for the name of the first scalar/vector/tensor DataArray. More...
 
std::string getEndian () const
 
std::string getFileExtension () const
 
Vtk::FormatTypes getFormat () const
 
Vtk::DataTypes getDatatype () const
 
auto comm () const
 

Protected Attributes

std::shared_ptr< DataCollectordataCollector_
 
Vtk::FormatTypes format_
 
Vtk::DataTypes datatype_
 
Vtk::DataTypes headertype_
 
Vtk::CompressorTypes compressor_ = Vtk::CompressorTypes::NONE
 
std::vector< VtkFunctionpointData_
 
std::vector< VtkFunctioncellData_
 
std::size_t const block_size = 1024*32
 
int compression_level = -1
 

Detailed Description

template<class GV, class DC>
class Dune::VtkWriterInterface< GV, DC >

Interface for file writers for the Vtk XML file formats.

Template Parameters
GVModel of Dune::GridView
DCModel of DataCollectorInterface

Member Typedef Documentation

◆ DataCollector

template<class GV , class DC >
using Dune::VtkWriterInterface< GV, DC >::DataCollector = DC

◆ GridView

template<class GV , class DC >
using Dune::VtkWriterInterface< GV, DC >::GridView = GV

◆ pos_type

template<class GV , class DC >
using Dune::VtkWriterInterface< GV, DC >::pos_type = typename std::ostream::pos_type
protected

◆ VtkFunction

template<class GV , class DC >
using Dune::VtkWriterInterface< GV, DC >::VtkFunction = Dune::Vtk::Function<GridView>
protected

Member Enumeration Documentation

◆ PositionTypes

template<class GV , class DC >
enum Dune::VtkWriterInterface::PositionTypes
protected
Enumerator
POINT_DATA 
CELL_DATA 

Constructor & Destructor Documentation

◆ VtkWriterInterface() [1/3]

template<class GV , class DC >
Dune::VtkWriterInterface< GV, DC >::VtkWriterInterface ( GridView const &  gridView,
Vtk::FormatTypes  format = Vtk::FormatTypes::BINARY,
Vtk::DataTypes  datatype = Vtk::DataTypes::FLOAT32,
Vtk::DataTypes  headertype = Vtk::DataTypes::UINT32 
)
inline

Constructor, passes the gridView to the DataCollector.

Creates a new VtkWriterInterface for the provided GridView. Initializes a DataCollector that is used to collect point coordinates, cell connectivity and data values.

This constructor assumes, that the DataCollector can be constructed from a single argument, the passed gridView.

Parameters
gridViewImplementation of Dune::GridView
formatFormat of the VTK file, either Vtk::FormatTypes::BINARY, Vtk::FormatTypes::ASCII, or Vtk::COMPRESSED
datatypeData type of a single component of the point coordinates [Vtk::DataTypes::FLOAT32]
headertypeInteger type used in binary data headers [Vtk::DataTypes::UINT32]

◆ VtkWriterInterface() [2/3]

template<class GV , class DC >
Dune::VtkWriterInterface< GV, DC >::VtkWriterInterface ( DataCollector dataCollector,
Vtk::FormatTypes  format = Vtk::FormatTypes::BINARY,
Vtk::DataTypes  datatype = Vtk::DataTypes::FLOAT32,
Vtk::DataTypes  headertype = Vtk::DataTypes::UINT32 
)
inline

Constructor, wraps the passed DataCollector in a non-destroying shared_ptr.

◆ VtkWriterInterface() [3/3]

template<class GV , class DC >
Dune::VtkWriterInterface< GV, DC >::VtkWriterInterface ( std::shared_ptr< DataCollector dataCollector,
Vtk::FormatTypes  format = Vtk::FormatTypes::BINARY,
Vtk::DataTypes  datatype = Vtk::DataTypes::FLOAT32,
Vtk::DataTypes  headertype = Vtk::DataTypes::UINT32 
)
inline

Constructor, stores the passed DataCollector.

Member Function Documentation

◆ addCellData()

template<class GV , class DC >
template<class Function , class... Args>
VtkWriterInterface& Dune::VtkWriterInterface< GV, DC >::addCellData ( Function &&  fct,
Args &&...  args 
)
inline

Attach cell data to the writer.

Attach a global function to the writer that will be evaluated at cell centers. The global function must be assignable to the function wrapper Vtk::Function. Additional argument for output datatype and number of components can be passed. See Vtk::Function Constructor for possible arguments.

Parameters
fctA GridFunction, LocalFunction, or Dune::VTKFunction
args...Additional arguments, like name, numComponents, dataType or Vtk::FieldInfo

◆ addPointData()

template<class GV , class DC >
template<class Function , class... Args>
VtkWriterInterface& Dune::VtkWriterInterface< GV, DC >::addPointData ( Function &&  fct,
Args &&...  args 
)
inline

Attach point data to the writer.

Attach a global function to the writer that will be evaluated at grid points (vertices and higher order points). The global function must be assignable to the function wrapper Vtk::Function. Additional argument for output datatype and number of components can be passed. See Vtk::Function Constructor for possible arguments.

Parameters
fctA GridFunction, LocalFunction, or Dune::VTKFunction
args...Additional arguments, like name, numComponents, dataType or Vtk::FieldInfo

◆ comm()

template<class GV , class DC >
auto Dune::VtkWriterInterface< GV, DC >::comm ( ) const
inlineprotected

◆ getDatatype()

template<class GV , class DC >
Vtk::DataTypes Dune::VtkWriterInterface< GV, DC >::getDatatype ( ) const
inlineprotected

◆ getEndian()

template<class GV , class DC >
std::string Dune::VtkWriterInterface< GV, DC >::getEndian ( ) const
inlineprotected

◆ getFileExtension()

template<class GV , class DC >
std::string Dune::VtkWriterInterface< GV, DC >::getFileExtension ( ) const
inlineprotected

◆ getFormat()

template<class GV , class DC >
Vtk::FormatTypes Dune::VtkWriterInterface< GV, DC >::getFormat ( ) const
inlineprotected

◆ getNames()

template<class GV , class DC >
std::string Dune::VtkWriterInterface< GV, DC >::getNames ( std::vector< VtkFunction > const &  data) const
protected

Return PointData/CellData attributes for the name of the first scalar/vector/tensor DataArray.

◆ setCompressor()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::setCompressor ( Vtk::CompressorTypes  compressor,
int  level = -1 
)
inline

Sets the compressor type used in binary data headers, Additionally a compression level can be passed with level = -1 means: default compression level. Level must be in [0-9]

◆ setDatatype()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::setDatatype ( Vtk::DataTypes  datatype)
inline

Sets the global datatype used for coordinates and other global float values.

◆ setFormat()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::setFormat ( Vtk::FormatTypes  format)
inline

◆ setHeadertype()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::setHeadertype ( Vtk::DataTypes  datatype)
inline

Sets the integer type used in binary data headers.

◆ write()

template<class GV , class DC >
std::string Dune::VtkWriterInterface< GV, DC >::write ( std::string const &  fn,
std::optional< std::string >  dir = {} 
) const
overridevirtual

Write the attached data to the file.

Parameters
fnFilename of the VTK file. May contain a directory and any file extension.
dirThe optional parameter specifies the directory of the partition files for parallel writes.
Returns
File name that is actually written.

Implements Dune::Vtk::FileWriter.

◆ writeAppended()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::writeAppended ( std::ofstream &  out,
std::vector< pos_type > const &  offsets 
) const
protected

◆ writeData()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::writeData ( std::ofstream &  out,
std::vector< pos_type > &  offsets,
VtkFunction const &  fct,
PositionTypes  type,
std::optional< std::size_t >  timestep = {} 
) const
protected

◆ writeDataAppended()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::writeDataAppended ( std::ofstream &  out,
std::vector< std::uint64_t > &  blocks 
) const
protected

◆ writeHeader()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::writeHeader ( std::ofstream &  out,
std::string const &  type 
) const
protected

◆ writePoints()

template<class GV , class DC >
void Dune::VtkWriterInterface< GV, DC >::writePoints ( std::ofstream &  out,
std::vector< pos_type > &  offsets,
std::optional< std::size_t >  timestep = {} 
) const
protected

◆ writeValuesAppended()

template<class GV , class DC >
template<class HeaderType , class FloatType >
std::uint64_t Dune::VtkWriterInterface< GV, DC >::writeValuesAppended ( std::ofstream &  out,
std::vector< FloatType > const &  values 
) const
protected

◆ writeValuesAscii()

template<class GV , class DC >
template<class T >
void Dune::VtkWriterInterface< GV, DC >::writeValuesAscii ( std::ofstream &  out,
std::vector< T > const &  values 
) const
protected

Member Data Documentation

◆ block_size

template<class GV , class DC >
std::size_t const Dune::VtkWriterInterface< GV, DC >::block_size = 1024*32
protected

◆ cellData_

template<class GV , class DC >
std::vector<VtkFunction> Dune::VtkWriterInterface< GV, DC >::cellData_
protected

◆ compression_level

template<class GV , class DC >
int Dune::VtkWriterInterface< GV, DC >::compression_level = -1
protected

◆ compressor_

template<class GV , class DC >
Vtk::CompressorTypes Dune::VtkWriterInterface< GV, DC >::compressor_ = Vtk::CompressorTypes::NONE
protected

◆ dataCollector_

template<class GV , class DC >
std::shared_ptr<DataCollector> Dune::VtkWriterInterface< GV, DC >::dataCollector_
protected

◆ datatype_

template<class GV , class DC >
Vtk::DataTypes Dune::VtkWriterInterface< GV, DC >::datatype_
protected

◆ format_

template<class GV , class DC >
Vtk::FormatTypes Dune::VtkWriterInterface< GV, DC >::format_
protected

◆ headertype_

template<class GV , class DC >
Vtk::DataTypes Dune::VtkWriterInterface< GV, DC >::headertype_
protected

◆ pointData_

template<class GV , class DC >
std::vector<VtkFunction> Dune::VtkWriterInterface< GV, DC >::pointData_
protected

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