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  Copyright (C) 2011-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef VTK_VECTOR_FUNCTION_HH
26 #define VTK_VECTOR_FUNCTION_HH
27 
29 
30 #include <dune/grid/io/file/vtk/function.hh>
31 #include <dune/istl/bvector.hh>
32 #include <dune/common/fvector.hh>
33 #include <dune/common/version.hh>
34 
35 #include <opm/common/Exceptions.hpp>
36 #include <opm/common/ErrorMacros.hpp>
37 
38 #include <string>
39 #include <limits>
40 #include <vector>
41 
42 namespace Ewoms {
43 
48 template <class GridView, class Mapper>
49 class VtkVectorFunction : public Dune::VTKFunction<GridView>
50 {
51  enum { dim = GridView::dimension };
52  typedef typename GridView::ctype ctype;
53  typedef typename GridView::template Codim<0>::Entity Element;
54 
55  typedef BaseOutputWriter::VectorBuffer VectorBuffer;
56 
57 public:
58  VtkVectorFunction(std::string name,
59  const GridView &gridView,
60  const Mapper &mapper,
61  const VectorBuffer &buf,
62  int codim)
63  : name_(name)
64  , gridView_(gridView)
65  , mapper_(mapper)
66  , buf_(buf)
67  , codim_(codim)
68  { assert(int(buf_.size()) == mapper_.size()); }
69 
70  virtual std::string name() const
71  { return name_; }
72 
73  virtual int ncomps() const
74  { return buf_[0].size(); }
75 
76  virtual double evaluate(int mycomp,
77  const Element &e,
78  const Dune::FieldVector<ctype, dim> &xi) const
79  {
80  int idx;
81  if (codim_ == 0) {
82  // cells. map element to the index
83 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
84  idx = mapper_.index(e);
85 #else
86  idx = mapper_.map(e);
87 #endif
88  }
89  else if (codim_ == dim) {
90  // find vertex which is closest to xi in local
91  // coordinates. This code is based on Dune::P1VTKFunction
92  double min = 1e100;
93  int imin = -1;
94  Dune::GeometryType gt = e.type();
95 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
96  int n = e.subEntities(dim);
97 #else
98  int n = e.template count<dim>();
99 #endif
100  for (int i = 0; i < n; ++i) {
101 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 3)
102  Dune::FieldVector<ctype, dim> local =
103  Dune::ReferenceElements<ctype, dim>::general(gt).position(i, dim);
104 #else
105  Dune::FieldVector<ctype, dim> local =
106  Dune::GenericReferenceElements<ctype, dim>::general(gt).position(i, dim);
107 #endif
108 
109  local -= xi;
110  if (local.infinity_norm() < min) {
111  min = local.infinity_norm();
112  imin = i;
113  }
114  }
115 
116  // map vertex to an index
117 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
118  idx = mapper_.subIndex(e, imin, codim_);
119 #else
120  idx = mapper_.map(e, imin, codim_);
121 #endif
122  }
123  else
124  OPM_THROW(std::logic_error, "Only element and vertex based vector "
125  " fields are supported so far.");
126 
127  double val = buf_[idx][mycomp];
128  if (std::abs(val) < std::numeric_limits<float>::min())
129  val = 0;
130 
131  return val;
132  }
133 
134 private:
135  const std::string name_;
136  const GridView gridView_;
137  const Mapper &mapper_;
138  const VectorBuffer &buf_;
139  int codim_;
140 };
141 
142 } // namespace Ewoms
143 
144 #endif
virtual int ncomps() const
Definition: vtkvectorfunction.hh:73
std::vector< Vector > VectorBuffer
Definition: baseoutputwriter.hh:48
Definition: baseauxiliarymodule.hh:35
The base class for all output writers.
virtual std::string name() const
Definition: vtkvectorfunction.hh:70
Provides a vector-valued function using Dune::FieldVectors as elements.
Definition: vtkvectorfunction.hh:49
virtual double evaluate(int mycomp, const Element &e, const Dune::FieldVector< ctype, dim > &xi) const
Definition: vtkvectorfunction.hh:76
VtkVectorFunction(std::string name, const GridView &gridView, const Mapper &mapper, const VectorBuffer &buf, int codim)
Definition: vtkvectorfunction.hh:58