25 #ifndef VTK_SCALAR_FUNCTION_HH
26 #define VTK_SCALAR_FUNCTION_HH
30 #include <dune/grid/io/file/vtk/function.hh>
31 #include <dune/istl/bvector.hh>
32 #include <dune/common/fvector.hh>
33 #include <dune/common/version.hh>
35 #include <opm/common/Exceptions.hpp>
36 #include <opm/common/ErrorMacros.hpp>
48 template <
class Gr
idView,
class Mapper>
51 enum { dim = GridView::dimension };
52 typedef typename GridView::ctype ctype;
53 typedef typename GridView::template Codim<0>::Entity Element;
59 const GridView &gridView,
61 const ScalarBuffer &buf,
68 { assert(
int(buf_.size()) == mapper_.size()); }
70 virtual std::string
name()
const
78 const Dune::FieldVector<ctype, dim> &xi)
const
83 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
84 idx = mapper_.index(e);
89 else if (codim_ == dim) {
94 Dune::GeometryType gt = e.type();
95 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
96 int n = e.subEntities(dim);
98 int n = e.template count<dim>();
100 for (
int i = 0; i < n; ++i) {
101 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 3)
102 Dune::FieldVector<ctype, dim> local =
103 Dune::ReferenceElements<ctype, dim>::general(gt).position(i, dim);
105 Dune::FieldVector<ctype, dim> local =
106 Dune::GenericReferenceElements<ctype, dim>::general(gt).position(i, dim);
110 if (local.infinity_norm() < min) {
111 min = local.infinity_norm();
117 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
118 idx = mapper_.subIndex(e, imin, codim_);
120 idx = mapper_.map(e, imin, codim_);
124 OPM_THROW(std::logic_error,
"Only element and vertex based vector "
125 " fields are supported so far.");
127 double val = buf_[idx];
128 if (std::abs(val) < std::numeric_limits<float>::min())
135 const std::string name_;
136 const GridView gridView_;
137 const Mapper &mapper_;
138 const ScalarBuffer &buf_;
std::vector< Scalar > ScalarBuffer
Definition: baseoutputwriter.hh:47
virtual double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtkscalarfunction.hh:76
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtkscalarfunction.hh:49
Definition: baseauxiliarymodule.hh:35
The base class for all output writers.
virtual std::string name() const
Definition: vtkscalarfunction.hh:70
VtkScalarFunction(std::string name, const GridView &gridView, const Mapper &mapper, const ScalarBuffer &buf, int codim)
Definition: vtkscalarfunction.hh:58
virtual int ncomps() const
Definition: vtkscalarfunction.hh:73