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