25 #ifndef VTK_TENSOR_FUNCTION_HH
26 #define VTK_TENSOR_FUNCTION_HH
30 #include <dune/grid/io/file/vtk/function.hh>
31 #include <dune/common/fvector.hh>
32 #include <dune/common/version.hh>
34 #include <opm/common/Exceptions.hpp>
35 #include <opm/common/ErrorMacros.hpp>
46 template <
class Gr
idView,
class Mapper>
49 enum { dim = GridView::dimension };
50 typedef typename GridView::ctype ctype;
51 typedef typename GridView::template Codim<0>::Entity Element;
57 const GridView &gridView,
59 const TensorBuffer &buf,
67 , matrixColumnIdx_(matrixColumnIdx)
68 { assert(
int(buf_.size()) == mapper_.size()); }
70 virtual std::string
name()
const
74 {
return buf_[0].M(); }
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,
125 "Only element and vertex based tensor fields are supported so far.");
128 int j = matrixColumnIdx_;
129 double val = buf_[idx][i][j];
130 if (std::abs(val) < std::numeric_limits<float>::min())
137 const std::string name_;
138 const GridView gridView_;
139 const Mapper &mapper_;
140 const TensorBuffer &buf_;
142 int matrixColumnIdx_;
virtual int ncomps() const
Definition: vtktensorfunction.hh:73
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:49
virtual double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtktensorfunction.hh:76
VtkTensorFunction(std::string name, const GridView &gridView, const Mapper &mapper, const TensorBuffer &buf, int codim, int matrixColumnIdx)
Definition: vtktensorfunction.hh:56
Definition: baseauxiliarymodule.hh:35
The base class for all output writers.
Provides a tensor-valued function using Dune::FieldMatrix objects as elements.
Definition: vtktensorfunction.hh:47
virtual std::string name() const
Definition: vtktensorfunction.hh:70