4 #ifndef DUNE_PDELAB_COMMON_VTKEXPORT_HH
5 #define DUNE_PDELAB_COMMON_VTKEXPORT_HH
12 #include<dune/common/shared_ptr.hh>
14 #include<dune/grid/io/file/vtk/vtkwriter.hh>
24 :
public Dune::VTKFunction<typename T::Traits::GridViewType>
26 typedef typename T::Traits::GridViewType::Grid::ctype DF;
27 enum {n=T::Traits::GridViewType::dimension};
28 typedef typename T::Traits::GridViewType::Grid::template Codim<0>::Entity Entity;
49 const std::vector<std::size_t> &remap_ =
51 : t(stackobject_to_shared_ptr(t_)), s(s_), remap(remap_)
69 const std::vector<std::size_t> &remap_ =
71 : t(t_), s(s_), remap(remap_)
79 virtual double evaluate (
int comp,
const Entity&
e,
const Dune::FieldVector<DF,n>& xi)
const override
81 typename T::Traits::DomainType x;
82 typename T::Traits::RangeType y;
84 for (
int i=0; i<n; i++)
87 return y[remap[comp]];
90 virtual std::string
name ()
const override
96 std::shared_ptr<const T> t;
98 std::vector<std::size_t> remap;
120 (
const std::shared_ptr<GF> &gf,
const std::string &name,
121 const std::vector<std::size_t> &remap =
123 {
return std::make_shared<VTKGridFunctionAdapter<GF> >(gf, name, remap); }
144 (
const GF &gf,
const std::string &name,
145 const std::vector<std::size_t> &remap =
147 {
return std::make_shared<VTKGridFunctionAdapter<GF> >(stackobject_to_shared_ptr(gf), name, remap); }
166 (
const std::shared_ptr<const GF> &gf,
const std::string &name,
167 const std::vector<std::size_t> &remap =
169 {
return std::make_shared<VTKGridFunctionAdapter<GF> >(gf, name, remap); }
175 template<
typename G,
typename FEM>
177 :
public Dune::VTKFunction<G>
179 typedef typename G::ctype DF;
180 enum {n=G::dimension};
181 typedef typename G::template Codim<0>::Entity Entity;
193 virtual double evaluate (
int comp,
const Entity&
e,
const Dune::FieldVector<DF,n>& xi)
const override
195 return fem.getOrder(
e);
198 virtual std::string
name ()
const override
const Entity & e
Definition: localfunctionspace.hh:121
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
std::vector< T > rangeVector(T begin, T passed_the_end, T increment=1)
Generate a vector with a range of numbers.
Definition: range.hh:22
std::shared_ptr< VTKGridFunctionAdapter< GF > > makeVTKGridFunctionAdapter(const std::shared_ptr< GF > &gf, const std::string &name, const std::vector< std::size_t > &remap=rangeVector(std::size_t(GF::Traits::dimRange)))
construct a VTKGridFunctionAdapter
Definition: vtkexport.hh:120
wrap a GridFunction so it can be used with the VTKWriter from dune-grid.
Definition: vtkexport.hh:25
virtual double evaluate(int comp, const Entity &e, const Dune::FieldVector< DF, n > &xi) const override
Definition: vtkexport.hh:79
virtual std::string name() const override
Definition: vtkexport.hh:90
VTKGridFunctionAdapter(const T &t_, std::string s_, const std::vector< std::size_t > &remap_=rangeVector(std::size_t(T::Traits::dimRange)))
construct a VTKGridFunctionAdapter
Definition: vtkexport.hh:48
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
Definition: vtkexport.hh:68
virtual int ncomps() const override
Definition: vtkexport.hh:74
Definition: vtkexport.hh:178
virtual std::string name() const override
Definition: vtkexport.hh:198
virtual double evaluate(int comp, const Entity &e, const Dune::FieldVector< DF, n > &xi) const override
Definition: vtkexport.hh:193
virtual int ncomps() const override
Definition: vtkexport.hh:188
VTKFiniteElementMapAdapter(const FEM &fem_, std::string s_)
Definition: vtkexport.hh:184