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
30#include <dune/common/fvector.hh>
31#include <dune/common/version.hh>
32
33#include <dune/grid/io/file/vtk/function.hh>
34
36
37#include <cstddef>
38#include <stdexcept>
39#include <string>
40#include <string_view>
41
42namespace Opm {
43
47template <class GridView, class Mapper>
48class VtkTensorFunction : public Dune::VTKFunction<GridView>
49{
50 enum { dim = GridView::dimension };
51 using ctype = typename GridView::ctype;
52 using Element = typename GridView::template Codim<0>::Entity;
53
54 using TensorBuffer = BaseOutputWriter::TensorBuffer;
55
56public:
57 VtkTensorFunction(std::string_view name,
58 const GridView& gridView,
59 const Mapper& mapper,
60 const TensorBuffer& buf,
61 unsigned codim,
62 unsigned matrixColumnIdx)
63 : name_(name)
64 , gridView_(gridView)
65 , mapper_(mapper)
66 , buf_(buf)
67 , codim_(codim)
68 , matrixColumnIdx_(matrixColumnIdx)
69 { assert(int(buf_.size()) == int(mapper_.size())); }
70
71 std::string name() const override
72 { return name_; }
73
74 int ncomps() const override
75 { return static_cast<int>(buf_[0].M()); }
76
77 double evaluate(int mycomp,
78 const Element& e,
79 const Dune::FieldVector<ctype, dim>& xi) const override
80 {
81 std::size_t idx;
82 if (codim_ == 0) {
83 // cells. map element to the index
84 idx = static_cast<std::size_t>(mapper_.index(e));
85 }
86 else if (codim_ == dim) {
87 // find vertex which is closest to xi in local
88 // coordinates. This code is based on Dune::P1VTKFunction
89 double min = 1e100;
90 int imin = -1;
91 Dune::GeometryType gt = e.type();
92 int n = static_cast<int>(e.subEntities(dim));
93 for (int i = 0; i < n; ++i) {
94 Dune::FieldVector<ctype, dim> local =
95 Dune::ReferenceElements<ctype, dim>::general(gt).position(i, dim);
96
97 local -= xi;
98 if (local.infinity_norm() < min) {
99 min = local.infinity_norm();
100 imin = i;
101 }
102 }
103
104 // map vertex to an index
105 idx = static_cast<std::size_t>(mapper_.subIndex(e, imin, codim_));
106 }
107 else {
108 throw std::logic_error("Only element and vertex based tensor fields are supported so far.");
109 }
110
111 const unsigned i = static_cast<unsigned>(mycomp);
112 const unsigned j = static_cast<unsigned>(matrixColumnIdx_);
113
114 return static_cast<double>(static_cast<float>(buf_[idx][i][j]));
115 }
116
117private:
118 const std::string name_;
119 const GridView gridView_;
120 const Mapper& mapper_;
121 const TensorBuffer& buf_;
122 const unsigned codim_;
123 const unsigned matrixColumnIdx_;
124};
125
126} // namespace Opm
127
128#endif
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:52
Provides a tensor-valued function using Dune::FieldMatrix objects as elements.
Definition: vtktensorfunction.hh:49
VtkTensorFunction(std::string_view name, const GridView &gridView, const Mapper &mapper, const TensorBuffer &buf, unsigned codim, unsigned matrixColumnIdx)
Definition: vtktensorfunction.hh:57
std::string name() const override
Definition: vtktensorfunction.hh:71
double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const override
Definition: vtktensorfunction.hh:77
int ncomps() const override
Definition: vtktensorfunction.hh:74
Definition: blackoilboundaryratevector.hh:39