27#ifndef VTK_TENSOR_FUNCTION_HH
28#define VTK_TENSOR_FUNCTION_HH
32#include <dune/grid/io/file/vtk/function.hh>
33#include <dune/common/fvector.hh>
34#include <dune/common/version.hh>
45template <
class Gr
idView,
class Mapper>
48 enum { dim = GridView::dimension };
49 using ctype =
typename GridView::ctype;
50 using Element =
typename GridView::template Codim<0>::Entity;
56 const GridView& gridView,
58 const TensorBuffer& buf,
60 unsigned matrixColumnIdx)
66 , matrixColumnIdx_(matrixColumnIdx)
67 { assert(
int(buf_.size()) ==
int(mapper_.size())); }
69 virtual std::string
name()
const
73 {
return static_cast<int>(buf_[0].M()); }
77 const Dune::FieldVector<ctype, dim>& xi)
const
82 idx =
static_cast<size_t>(mapper_.index(e));
84 else if (codim_ == dim) {
89 Dune::GeometryType gt = e.type();
90 int n =
static_cast<int>(e.subEntities(dim));
91 for (
int i = 0; i < n; ++i) {
92 Dune::FieldVector<ctype, dim> local =
93 Dune::ReferenceElements<ctype, dim>::general(gt).position(i, dim);
96 if (local.infinity_norm() < min) {
97 min = local.infinity_norm();
103 idx =
static_cast<size_t>(mapper_.subIndex(e, imin, codim_));
106 throw std::logic_error(
"Only element and vertex based tensor fields are supported so far.");
108 unsigned i =
static_cast<unsigned>(mycomp);
109 unsigned j =
static_cast<unsigned>(matrixColumnIdx_);
111 return static_cast<double>(
static_cast<float>(buf_[idx][i][j]));
115 const std::string name_;
116 const GridView gridView_;
117 const Mapper& mapper_;
118 const TensorBuffer& buf_;
120 unsigned matrixColumnIdx_;
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:51
Provides a tensor-valued function using Dune::FieldMatrix objects as elements.
Definition: vtktensorfunction.hh:47
virtual int ncomps() const
Definition: vtktensorfunction.hh:72
virtual std::string name() const
Definition: vtktensorfunction.hh:69
VtkTensorFunction(std::string name, const GridView &gridView, const Mapper &mapper, const TensorBuffer &buf, unsigned codim, unsigned matrixColumnIdx)
Definition: vtktensorfunction.hh:55
virtual double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtktensorfunction.hh:75
Definition: blackoilboundaryratevector.hh:37