vtktensorfunction.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef VTK_TENSOR_FUNCTION_HH
28#define VTK_TENSOR_FUNCTION_HH
29
31
32#include <dune/grid/io/file/vtk/function.hh>
33#include <dune/common/fvector.hh>
34#include <dune/common/version.hh>
35
36#include <string>
37#include <limits>
38#include <vector>
39
40namespace Opm {
41
45template <class GridView, class Mapper>
46class VtkTensorFunction : public Dune::VTKFunction<GridView>
47{
48 enum { dim = GridView::dimension };
49 using ctype = typename GridView::ctype;
50 using Element = typename GridView::template Codim<0>::Entity;
51
52 using TensorBuffer = BaseOutputWriter::TensorBuffer;
53
54public:
56 const GridView& gridView,
57 const Mapper& mapper,
58 const TensorBuffer& buf,
59 unsigned codim,
60 unsigned matrixColumnIdx)
61 : name_(name)
62 , gridView_(gridView)
63 , mapper_(mapper)
64 , buf_(buf)
65 , codim_(codim)
66 , matrixColumnIdx_(matrixColumnIdx)
67 { assert(int(buf_.size()) == int(mapper_.size())); }
68
69 virtual std::string name() const
70 { return name_; }
71
72 virtual int ncomps() const
73 { return static_cast<int>(buf_[0].M()); }
74
75 virtual double evaluate(int mycomp,
76 const Element& e,
77 const Dune::FieldVector<ctype, dim>& xi) const
78 {
79 size_t idx;
80 if (codim_ == 0) {
81 // cells. map element to the index
82 idx = static_cast<size_t>(mapper_.index(e));
83 }
84 else if (codim_ == dim) {
85 // find vertex which is closest to xi in local
86 // coordinates. This code is based on Dune::P1VTKFunction
87 double min = 1e100;
88 int imin = -1;
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);
94
95 local -= xi;
96 if (local.infinity_norm() < min) {
97 min = local.infinity_norm();
98 imin = i;
99 }
100 }
101
102 // map vertex to an index
103 idx = static_cast<size_t>(mapper_.subIndex(e, imin, codim_));
104 }
105 else
106 throw std::logic_error("Only element and vertex based tensor fields are supported so far.");
107
108 unsigned i = static_cast<unsigned>(mycomp);
109 unsigned j = static_cast<unsigned>(matrixColumnIdx_);
110
111 return static_cast<double>(static_cast<float>(buf_[idx][i][j]));
112 }
113
114private:
115 const std::string name_;
116 const GridView gridView_;
117 const Mapper& mapper_;
118 const TensorBuffer& buf_;
119 unsigned codim_;
120 unsigned matrixColumnIdx_;
121};
122
123} // namespace Opm
124
125#endif
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