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