25#ifndef VTK_NESTED_FUNCTION_HH
26#define VTK_NESTED_FUNCTION_HH
30#include <dune/common/fvector.hh>
31#include <dune/istl/bvector.hh>
32#include <dune/grid/io/file/vtk/function.hh>
41template <
class Gr
idView,
class Mapper,
class Buffer>
44 enum { dim = GridView::dimension };
45 using ctype =
typename GridView::ctype;
46 using Element =
typename GridView::template Codim<0>::Entity;
47 using ReferenceElements = Dune::ReferenceElements<ctype, dim>;
51 const GridView &gridView,
63 assert(std::size_t(buf_.size()) == std::size_t(mapper_.size()));
66 virtual std::string
name ()
const
73 const Element &element,
74 const Dune::FieldVector< ctype, dim > &xi)
const
79 idx = mapper_.index(element);
81 else if (codim_ == dim) {
86 Dune::GeometryType geomType = element.type();
87 int n = element.subEntities(dim);
89 for (
int i=0; i < n; ++i)
91 Dune::FieldVector<ctype,dim> local =
92 ReferenceElements::general(geomType).position(i,dim);
94 if (local.infinity_norm()<min)
96 min = local.infinity_norm();
102 idx = mapper_.subIndex(element, imin, codim_);
105 DUNE_THROW(Dune::InvalidStateException,
106 "Only element and vertex based vector "
107 " fields are supported so far.");
109 double val = buf_[idx][mycomp];
111 if (abs(val) < std::numeric_limits<float>::min())
118 const std::string name_;
119 const GridView gridView_;
120 const Mapper &mapper_;
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtknestedfunction.hh:43
virtual int ncomps() const
Definition: vtknestedfunction.hh:69
virtual double evaluate(int mycomp, const Element &element, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtknestedfunction.hh:72
virtual std::string name() const
Definition: vtknestedfunction.hh:66
VtkNestedFunction(std::string name, const GridView &gridView, const Mapper &mapper, const Buffer &buf, int codim, int numComp)
Definition: vtknestedfunction.hh:50