DUNE PDELab (git)

Dune::PDELab::VTKGridFunctionAdapter< T > Class Template Referenceabstract

wrap a GridFunction so it can be used with the VTKWriter from dune-grid. More...

#include <dune/pdelab/common/vtkexport.hh>

Public Member Functions

 VTKGridFunctionAdapter (const T &t_, std::string s_, const std::vector< std::size_t > &remap_=rangeVector(std::size_t(T::Traits::dimRange)))
 construct a VTKGridFunctionAdapter More...
 
 VTKGridFunctionAdapter (const std::shared_ptr< const T > &t_, std::string s_, const std::vector< std::size_t > &remap_=rangeVector(std::size_t(T::Traits::dimRange)))
 construct a VTKGridFunctionAdapter More...
 
virtual int ncomps () const override
 
virtual std::string name () const override
 get name
 
virtual double evaluate (int comp, const Entity &e, const Dune::FieldVector< ctype, dim > &xi) const=0
 evaluate single component comp in the entity e at local coordinates xi More...
 
virtual VTK::Precision precision () const
 get output precision for the field
 

Detailed Description

template<typename T>
class Dune::PDELab::VTKGridFunctionAdapter< T >

wrap a GridFunction so it can be used with the VTKWriter from dune-grid.

Examples
recipe-communication.cc.

Constructor & Destructor Documentation

◆ VTKGridFunctionAdapter() [1/2]

template<typename T >
Dune::PDELab::VTKGridFunctionAdapter< T >::VTKGridFunctionAdapter ( const T &  t_,
std::string  s_,
const std::vector< std::size_t > &  remap_ = rangeVector(std::size_t(T::Traits::dimRange)) 
)
inline

construct a VTKGridFunctionAdapter

Parameters
t_GridFunction object to wrap. A reference to the grid function object is stored internally and the constructed object becomes invalid as soon as that reference becomes invalid.
s_Name of the field as returned by name().
remap_How components are remapped between PDELab and VTK. The default value yields the identity map with entries (0,1,...,T::Traits::dimRange-1).

The resulting VTKFunction will have remap_.size() components. None of the elements of remap_ should be greater or equal to T::Traits::dimRange. When component c is requested by the VTKWriter, it will be mapped to the component remap_[c] of the GridFunction.

◆ VTKGridFunctionAdapter() [2/2]

template<typename T >
Dune::PDELab::VTKGridFunctionAdapter< T >::VTKGridFunctionAdapter ( const std::shared_ptr< const T > &  t_,
std::string  s_,
const std::vector< std::size_t > &  remap_ = rangeVector(std::size_t(T::Traits::dimRange)) 
)
inline

construct a VTKGridFunctionAdapter

Parameters
t_Shared pointer to a GridFunction object to wrap.
s_Name of the field as returned by name().
remap_How components are remapped between PDELab and VTK. The default value yields the identity map with entries (0,1,...,T::Traits::dimRange-1).

The resulting VTKFunction will have remap_.size() components. None of the elements of remap_ should be greater or equal to T::Traits::dimRange. When component c is requested by the VTKWriter, it will be mapped to the component remap_[c] of the GridFunction.

Member Function Documentation

◆ evaluate()

virtual double Dune::VTKFunction< T::Traits::GridViewType >::evaluate ( int  comp,
const Entity &  e,
const Dune::FieldVector< ctype, dim > &  xi 
) const
pure virtualinherited

evaluate single component comp in the entity e at local coordinates xi

Evaluate the function in an entity at local coordinates.

Parameters
[in]compnumber of component to be evaluated
[in]ereference to grid entity of codimension 0
[in]xipoint in local coordinates of the reference element of e
Returns
value of the component

◆ ncomps()

template<typename T >
virtual int Dune::PDELab::VTKGridFunctionAdapter< T >::ncomps ( ) const
inlineoverridevirtual

return number of components (1 for scalar valued functions, 3 for vector valued function in 3D etc.)

Implements Dune::VTKFunction< T::Traits::GridViewType >.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 23, 23:29, 2024)