template<class GridView, Dune::VTK::OutputType OutputValue = Dune::VTK::ascii>
class Dumux::VtkMultiWriter< GridView, OutputValue >
Simplifies writing multi-file VTK datasets.
This class automatically keeps the meta file up to date and simplifies writing datasets consisting of multiple files. (i.e. multiple time steps or grid refinements within a time step.)
- Todo:
- This class can most likely be replaced by Dune::VTKSequenceWriter
|
| VtkMultiWriter (const GridView &gridView, const std::string &simName="", std::string multiFileName="") |
|
| ~VtkMultiWriter () |
|
int | curWriterNum () const |
| Returns the number of the current VTK file. More...
|
|
void | gridChanged () |
| Updates the internal data structures after mesh refinement. More...
|
|
void | beginWrite (double t) |
| Called when ever a new time step or a new grid must be written. More...
|
|
template<class Scalar , int nComp> |
Dune::BlockVector< Dune::FieldVector< Scalar, nComp > > * | allocateManagedBuffer (int nEntities) |
| Allocate a managed buffer for a scalar field. More...
|
|
template<class Scalar > |
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > * | allocateManagedBuffer (int nEntities) |
|
Dune::BlockVector< Dune::FieldVector< double, 1 > > * | allocateManagedBuffer (int nEntities) |
|
template<class DataBuffer > |
void | attachVertexData (DataBuffer &buf, std::string name, int nComps=1) |
| Add a finished vertex centered vector field to the output. More...
|
|
template<class DataBuffer > |
void | attachCellData (DataBuffer &buf, std::string name, int nComps=1) |
| Add a cell centered quantity to the output. More...
|
|
template<class DataBuffer > |
void | attachDofData (DataBuffer &buf, std::string name, bool isVertexData, int nComps=1) |
| Add data associated with degrees of freedom to the output. More...
|
|
void | endWrite (bool onlyDiscard=false) |
| Finalizes the current writer. More...
|
|
template<class Restarter > |
void | serialize (Restarter &res) |
| Write the multi-writer's state to a restart file. More...
|
|
template<class Restarter > |
void | deserialize (Restarter &res) |
| Read the multi-writer's state from a restart file. More...
|
|
template<class GridView , Dune::VTK::OutputType OutputValue = Dune::VTK::ascii>
template<class Scalar , int nComp>
Dune::BlockVector< Dune::FieldVector< Scalar, nComp > > * Dumux::VtkMultiWriter< GridView, OutputValue >::allocateManagedBuffer |
( |
int |
nEntities | ) |
|
|
inline |
Allocate a managed buffer for a scalar field.
The buffer will be deleted automatically after the data has been written by to disk.
template<class GridView , Dune::VTK::OutputType OutputValue = Dune::VTK::ascii>
template<class DataBuffer >
void Dumux::VtkMultiWriter< GridView, OutputValue >::attachCellData |
( |
DataBuffer & |
buf, |
|
|
std::string |
name, |
|
|
int |
nComps = 1 |
|
) |
| |
|
inline |
Add a cell centered quantity to the output.
If the buffer is managed by the VtkMultiWriter, it must have been created using createField() and may not be used by anywhere after calling this method. After the data is written to disk, it will be deleted automatically.
If the buffer is not managed by the MultiWriter, the buffer must exist at least until the call to endWrite() finishes.
In both cases, modifying the buffer between the call to this method and endWrite() results in undefined behaviour.
template<class GridView , Dune::VTK::OutputType OutputValue = Dune::VTK::ascii>
template<class DataBuffer >
void Dumux::VtkMultiWriter< GridView, OutputValue >::attachDofData |
( |
DataBuffer & |
buf, |
|
|
std::string |
name, |
|
|
bool |
isVertexData, |
|
|
int |
nComps = 1 |
|
) |
| |
|
inline |
Add data associated with degrees of freedom to the output.
This is a wrapper for the functions attachVertexData and attachCelldata. Depending on the value of isVertexData, the appropriate function is selected.
template<class GridView , Dune::VTK::OutputType OutputValue = Dune::VTK::ascii>
template<class DataBuffer >
void Dumux::VtkMultiWriter< GridView, OutputValue >::attachVertexData |
( |
DataBuffer & |
buf, |
|
|
std::string |
name, |
|
|
int |
nComps = 1 |
|
) |
| |
|
inline |
Add a finished vertex centered vector field to the output.
Add a vertex-centered quantity to the output.
If the buffer is managed by the VtkMultiWriter, it must have been created using createField() and may not be used by anywhere after calling this method. After the data is written to disk, it will be deleted automatically.
If the buffer is not managed by the MultiWriter, the buffer must exist at least until the call to endWrite() finishes.
In both cases, modifying the buffer between the call to this method and endWrite() results in undefined behaviour.
template<class GridView , Dune::VTK::OutputType OutputValue = Dune::VTK::ascii>
Finalizes the current writer.
This means that everything will be written to disk, except if the onlyDiscard argument is true. In this case only all managed buffers are deleted, but no output is written.