opm-simulators
Loading...
Searching...
No Matches
Opm::VtkMultiWriter< GridView, vtkFormat > Class Template Reference

Simplifies writing multi-file VTK datasets. More...

#include <vtkmultiwriter.hh>

Inheritance diagram for Opm::VtkMultiWriter< GridView, vtkFormat >:
Opm::BaseOutputWriter

Public Types

using Scalar = BaseOutputWriter::Scalar
using Tensor = BaseOutputWriter::Tensor
using Vector = BaseOutputWriter::Vector
using ScalarBuffer = BaseOutputWriter::ScalarBuffer
using TensorBuffer = BaseOutputWriter::TensorBuffer
using VectorBuffer = BaseOutputWriter::VectorBuffer
using VtkWriter = Dune::VTKWriter<GridView>
Public Types inherited from Opm::BaseOutputWriter
using Scalar = double
using Tensor = Dune::DynamicMatrix<double>
using Vector = Dune::DynamicVector<double>
using ScalarBuffer = std::vector<Scalar>
using TensorBuffer = std::vector<Tensor>
using VectorBuffer = std::vector<Vector>

Public Member Functions

 VtkMultiWriter (bool asyncWriting, const GridView &gridView, const std::string &outputDir, const std::string &simName="", const std::string &multiFileName="")
int curWriterNum () const
 Returns the number of the current VTK file.
void gridChanged ()
 Updates the internal data structures after mesh refinement.
void beginWrite (double t) override
 Called whenever a new time step must be written.
ScalarBuffer * allocateManagedScalarBuffer (std::size_t numEntities)
 Allocate a managed buffer for a scalar field.
VectorBuffer * allocateManagedVectorBuffer (std::size_t numOuter, std::size_t numInner)
 Allocate a managed buffer for a vector field.
void attachScalarVertexData (ScalarBuffer &buf, std::string_view name) override
 Add a finished vertex centered vector field to the output.
void attachScalarElementData (ScalarBuffer &buf, std::string_view name) override
 Add a element centered quantity to the output.
void attachVectorVertexData (VectorBuffer &buf, std::string_view name) override
 Add a finished vertex centered vector field to the output.
void attachTensorVertexData (TensorBuffer &buf, std::string_view name) override
 Add a finished vertex-centered tensor field to the output.
void attachVectorElementData (VectorBuffer &buf, std::string_view name) override
 Add a element centered quantity to the output.
void attachTensorElementData (TensorBuffer &buf, std::string_view name) override
 Add a finished element-centered tensor field to the output.
void endWrite (bool onlyDiscard) override
 Finalizes the current writer.
template<class Restarter>
void serialize (Restarter &res)
 Write the multi-writer's state to a restart file.
template<class Restarter>
void deserialize (Restarter &res)
 Read the multi-writer's state from a restart file.

Detailed Description

template<class GridView, Dune::VTK::OutputType vtkFormat>
class Opm::VtkMultiWriter< GridView, vtkFormat >

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.)

Member Function Documentation

◆ allocateManagedScalarBuffer()

template<class GridView, Dune::VTK::OutputType vtkFormat>
ScalarBuffer * Opm::VtkMultiWriter< GridView, vtkFormat >::allocateManagedScalarBuffer ( std::size_t numEntities)
inline

Allocate a managed buffer for a scalar field.

The buffer will be deleted automatically after the data has been written by to disk.

◆ allocateManagedVectorBuffer()

template<class GridView, Dune::VTK::OutputType vtkFormat>
VectorBuffer * Opm::VtkMultiWriter< GridView, vtkFormat >::allocateManagedVectorBuffer ( std::size_t numOuter,
std::size_t numInner )
inline

Allocate a managed buffer for a vector field.

The buffer will be deleted automatically after the data has been written by to disk.

◆ attachScalarElementData()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::attachScalarElementData ( ScalarBuffer & buf,
std::string_view name )
inlineoverridevirtual

Add a element 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.

Implements Opm::BaseOutputWriter.

◆ attachScalarVertexData()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::attachScalarVertexData ( ScalarBuffer & buf,
std::string_view name )
inlineoverridevirtual

Add a finished vertex centered vector field to the output.

If the buffer is managed by the VtkMultiWriter, it must have been created using allocateManagedBuffer() and may not be used 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 behavior.

Implements Opm::BaseOutputWriter.

◆ attachTensorElementData()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::attachTensorElementData ( TensorBuffer & buf,
std::string_view name )
inlineoverridevirtual

Add a finished element-centered tensor field to the output.

Implements Opm::BaseOutputWriter.

◆ attachTensorVertexData()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::attachTensorVertexData ( TensorBuffer & buf,
std::string_view name )
inlineoverridevirtual

Add a finished vertex-centered tensor field to the output.

Implements Opm::BaseOutputWriter.

◆ attachVectorElementData()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::attachVectorElementData ( VectorBuffer & buf,
std::string_view name )
inlineoverridevirtual

Add a element 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.

Implements Opm::BaseOutputWriter.

◆ attachVectorVertexData()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::attachVectorVertexData ( VectorBuffer & buf,
std::string_view name )
inlineoverridevirtual

Add a finished vertex centered vector field to the output.

If the buffer is managed by the VtkMultiWriter, it must have been created using allocateManagedBuffer() and may not be used 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 behavior.

Implements Opm::BaseOutputWriter.

◆ beginWrite()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::beginWrite ( double t)
inlineoverridevirtual

Called whenever a new time step must be written.

Implements Opm::BaseOutputWriter.

◆ endWrite()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::endWrite ( bool onlyDiscard)
inlineoverridevirtual

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.

Implements Opm::BaseOutputWriter.

◆ gridChanged()

template<class GridView, Dune::VTK::OutputType vtkFormat>
void Opm::VtkMultiWriter< GridView, vtkFormat >::gridChanged ( )
inline

Updates the internal data structures after mesh refinement.

If the grid changes between two calls of beginWrite(), this method must be called before the second beginWrite()!


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