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