dune-vtk  0.2
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Protected Member Functions | List of all members
Dune::VtkReader< Grid, GC, FieldType > Class Template Reference

File-Reader for Vtk unstructured .vtu files. More...

#include <dune/vtk/vtkreader.hh>

Inheritance diagram for Dune::VtkReader< Grid, GC, FieldType >:
Inheritance graph

Public Types

using GridCreator = GC
 
using PointGridFunction = GridFunction< Vtk::PointContext >
 GridFunction representing the data stored on the points in the file. More...
 
using CellGridFunction = GridFunction< Vtk::CellContext >
 GridFunction representing the data stored on the cells in the file. More...
 

Public Member Functions

template<class... Args, std::enable_if_t< std::is_constructible< GridCreator, Args... >::value, int > = 0>
 VtkReader (Args &&... args)
 Constructor. Creates a new GridCreator with the passed factory. More...
 
 VtkReader (GridCreator &creator)
 Constructor. Stores the references in a non-destroying shared_ptr. More...
 
 VtkReader (std::shared_ptr< GridCreator > creator)
 Constructor. Stores the shared_ptr. More...
 
void read (std::string const &filename, bool fillCreator=true)
 Read the grid from file with filename into the GridCreator. More...
 
GridCreatorgridCreator ()
 Obtains the creator of the reader. More...
 
std::unique_ptr< Grid > createGrid () const
 
GridFunction< Vtk::PointContextgetPointData (std::string const &name) const
 
std::vector< DataArrayAttributes > getPointDataAttributes () const
 
GridFunction< Vtk::CellContextgetCellData (std::string const &name) const
 
std::vector< DataArrayAttributes > getCellDataAttributes () const
 
std::vector< std::string > const & pieces () const
 Return the filenames of parallel pieces. More...
 
template<class FloatType , class HeaderType >
void readCellDataAppended (MetaType< FloatType >, MetaType< HeaderType >, std::ifstream &input, std::string id)
 
template<class FloatType , class HeaderType >
void readPointDataAppended (MetaType< FloatType >, MetaType< HeaderType >, std::ifstream &input, std::string id)
 
template<class FloatType , class HeaderType >
void readPointsAppended (MetaType< FloatType >, MetaType< HeaderType >, std::ifstream &input)
 
template<class HeaderType >
void readCellsAppended (MetaType< HeaderType >, std::ifstream &input)
 
void readSerialFileFromStream (std::ifstream &input, bool create=true)
 Read the grid from an input stream, referring to a .vtu file, into the GridFactory factory_. More...
 
void readParallelFileFromStream (std::ifstream &input, int rank, int size, bool create=true)
 Read the grid from and input stream, referring to a .pvtu file, into the GridFactory factory_. More...
 
void fillGridCreator (bool insertPieces=true)
 

Static Public Member Functions

static std::unique_ptr< Grid > createGridFromFile (const std::string &filename, Args &&... args)
 
static void fillFactory (GridFactory< Grid > &factory, const std::string &filename, Args &&... args)
 

Static Protected Member Functions

static std::unique_ptr< Grid > createGridFromFileImpl (const std::string &filename, Args &&... args)
 
static void fillFactoryImpl (GridFactory< Grid > &, const std::string &, Args &&...)
 

Detailed Description

template<class Grid, class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
class Dune::VtkReader< Grid, GC, FieldType >

File-Reader for Vtk unstructured .vtu files.

Reads .vtu files and constructs a grid from the cells stored in the file Additionally, stored data can be read.

NOTE: Assumption on the file structure: Each XML tag must be on a separate line.

Template Parameters
GridThe type of the grid to construct.
GCGridCreator policy type to control what to pass to a grid factory with data given from the file. [ContinuousGridCreator]
FieldTypeType of the components of the data to extract from the file [default: double]

Member Typedef Documentation

◆ CellGridFunction

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
using Dune::VtkReader< Grid, GC, FieldType >::CellGridFunction = GridFunction<Vtk::CellContext>

GridFunction representing the data stored on the cells in the file.

◆ GridCreator

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
using Dune::VtkReader< Grid, GC, FieldType >::GridCreator = GC

◆ PointGridFunction

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
using Dune::VtkReader< Grid, GC, FieldType >::PointGridFunction = GridFunction<Vtk::PointContext>

GridFunction representing the data stored on the points in the file.

Constructor & Destructor Documentation

◆ VtkReader() [1/3]

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
template<class... Args, std::enable_if_t< std::is_constructible< GridCreator, Args... >::value, int > = 0>
Dune::VtkReader< Grid, GC, FieldType >::VtkReader ( Args &&...  args)
inlineexplicit

Constructor. Creates a new GridCreator with the passed factory.

Parameters
args...Either pass a GridFactory by reference or shared_ptr, or a list of arguments passed to the constructor of a Dune::GridFactory (typically and empty parameter list). See the constructor of GridCreatorInterface and the GridCreator passed to this reader.

◆ VtkReader() [2/3]

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
Dune::VtkReader< Grid, GC, FieldType >::VtkReader ( GridCreator creator)
inlineexplicit

Constructor. Stores the references in a non-destroying shared_ptr.

◆ VtkReader() [3/3]

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
Dune::VtkReader< Grid, GC, FieldType >::VtkReader ( std::shared_ptr< GridCreator creator)
inlineexplicit

Constructor. Stores the shared_ptr.

Member Function Documentation

◆ createGrid()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
std::unique_ptr<Grid> Dune::VtkReader< Grid, GC, FieldType >::createGrid ( ) const
inline

Construct the actual grid using the GridCreator [[expects: read_ == true]]

◆ createGridFromFile()

static std::unique_ptr<Grid> Dune::Vtk::FileReader< Grid, VtkReader< Grid, Vtk::ContinuousGridCreator< Grid > > >::createGridFromFile ( const std::string &  filename,
Args &&...  args 
)
inlinestaticinherited

Reads the grid from a file with filename and returns a unique_ptr to the created grid. Redirects to concrete implementation of derivated class.

◆ createGridFromFileImpl()

static std::unique_ptr<Grid> Dune::Vtk::FileReader< Grid, VtkReader< Grid, Vtk::ContinuousGridCreator< Grid > > >::createGridFromFileImpl ( const std::string &  filename,
Args &&...  args 
)
inlinestaticprotectedinherited

◆ fillFactory()

static void Dune::Vtk::FileReader< Grid, VtkReader< Grid, Vtk::ContinuousGridCreator< Grid > > >::fillFactory ( GridFactory< Grid > &  factory,
const std::string &  filename,
Args &&...  args 
)
inlinestaticinherited

Reads the grid from a file with filename into a grid-factory. Redirects to concrete implementation of derivated class.

◆ fillFactoryImpl()

static void Dune::Vtk::FileReader< Grid, VtkReader< Grid, Vtk::ContinuousGridCreator< Grid > > >::fillFactoryImpl ( GridFactory< Grid > &  ,
const std::string &  ,
Args &&  ... 
)
inlinestaticprotectedinherited

◆ fillGridCreator()

template<class Grid , class Creator , class Field >
void Dune::VtkReader< Grid, Creator, Field >::fillGridCreator ( bool  insertPieces = true)

Insert all internal data to the GridCreator NOTE: requires an aforegoing call to read()

◆ getCellData()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
GridFunction<Vtk::CellContext> Dune::VtkReader< Grid, GC, FieldType >::getCellData ( std::string const &  name) const
inline

Construct a grid-function representing the cell-data with the given name [[expects: read_ == true]]

◆ getCellDataAttributes()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
std::vector<DataArrayAttributes> Dune::VtkReader< Grid, GC, FieldType >::getCellDataAttributes ( ) const
inline

Return a vector of DataArrayAttributes for all CELL_DATA blocks [[expects: read_ == true]]

◆ getPointData()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
GridFunction<Vtk::PointContext> Dune::VtkReader< Grid, GC, FieldType >::getPointData ( std::string const &  name) const
inline

Construct a grid-function representing the point-data with the given name [[expects: read_ == true]]

◆ getPointDataAttributes()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
std::vector<DataArrayAttributes> Dune::VtkReader< Grid, GC, FieldType >::getPointDataAttributes ( ) const
inline

Return a vector of DataArrayAttributes for all POINT_DATA blocks [[expects: read_ == true]]

◆ gridCreator()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
GridCreator& Dune::VtkReader< Grid, GC, FieldType >::gridCreator ( )
inline

Obtains the creator of the reader.

◆ pieces()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
std::vector<std::string> const& Dune::VtkReader< Grid, GC, FieldType >::pieces ( ) const
inline

Return the filenames of parallel pieces.

◆ read()

template<class Grid , class Creator , class Field >
void Dune::VtkReader< Grid, Creator, Field >::read ( std::string const &  filename,
bool  fillCreator = true 
)

Read the grid from file with filename into the GridCreator.

This function fills internal data containers representing the information from the passed file.

Parameters
filenameThe name of the input file
fillCreatorIf false, only fill internal data structures, if true, pass the internal data to the GridCreator. [true]

◆ readCellDataAppended()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
template<class FloatType , class HeaderType >
void Dune::VtkReader< Grid, GC, FieldType >::readCellDataAppended ( MetaType< FloatType >  ,
MetaType< HeaderType >  ,
std::ifstream &  input,
std::string  id 
)

◆ readCellsAppended()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
template<class HeaderType >
void Dune::VtkReader< Grid, GC, FieldType >::readCellsAppended ( MetaType< HeaderType >  ,
std::ifstream &  input 
)

◆ readParallelFileFromStream()

template<class Grid , class Creator , class Field >
void Dune::VtkReader< Grid, Creator, Field >::readParallelFileFromStream ( std::ifstream &  input,
int  rank,
int  size,
bool  create = true 
)

Read the grid from and input stream, referring to a .pvtu file, into the GridFactory factory_.

Parameters
inputA STL input stream to read the VTK file from.
createIf false, only fill internal data structures, if true, also create the grid. [true]

◆ readPointDataAppended()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
template<class FloatType , class HeaderType >
void Dune::VtkReader< Grid, GC, FieldType >::readPointDataAppended ( MetaType< FloatType >  ,
MetaType< HeaderType >  ,
std::ifstream &  input,
std::string  id 
)

◆ readPointsAppended()

template<class Grid , class GC = Vtk::ContinuousGridCreator<Grid>, class FieldType = double>
template<class FloatType , class HeaderType >
void Dune::VtkReader< Grid, GC, FieldType >::readPointsAppended ( MetaType< FloatType >  ,
MetaType< HeaderType >  ,
std::ifstream &  input 
)

◆ readSerialFileFromStream()

template<class Grid , class Creator , class Field >
void Dune::VtkReader< Grid, Creator, Field >::readSerialFileFromStream ( std::ifstream &  input,
bool  create = true 
)

Read the grid from an input stream, referring to a .vtu file, into the GridFactory factory_.

Advanced read methods

Parameters
inputA STL input stream to read the VTK file from.
createIf false, only fill internal data structures, if true, also create the grid. [true]

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