NonuniformTableLinear.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2009, 2010, 2011, 2012 SINTEF ICT, Applied Mathematics.
3  Copyright 2009, 2010, 2011, 2012 Statoil ASA.
4 
5  This file is part of the Open Porous Media project (OPM).
6 
7  OPM is free software: you can redistribute it and/or modify
8  it under the terms of the GNU General Public License as published by
9  the Free Software Foundation, either version 3 of the License, or
10  (at your option) any later version.
11 
12  OPM is distributed in the hope that it will be useful,
13  but WITHOUT ANY WARRANTY; without even the implied warranty of
14  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  GNU General Public License for more details.
16 
17  You should have received a copy of the GNU General Public License
18  along with OPM. If not, see <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef OPM_NONUNIFORMTABLELINEAR_HEADER_INCLUDED
22 #define OPM_NONUNIFORMTABLELINEAR_HEADER_INCLUDED
23 
24 #include <cmath>
25 #include <exception>
26 #include <vector>
27 #include <utility>
28 
29 #include <opm/common/ErrorMacros.hpp>
31 
32 namespace Opm
33 {
34 
35 
41  template<typename T>
43  {
44  public:
47 
51  NonuniformTableLinear(const std::vector<double>& x_values,
52  const std::vector<T>& y_values);
53 
56  std::pair<double, double> domain();
57 
60  void rescaleDomain(std::pair<double, double> new_domain);
61 
65  double operator()(const double x) const;
66 
70  double derivative(const double x) const;
71 
75  double inverse(const double y) const;
76 
80  bool operator==(const NonuniformTableLinear& other) const;
81 
82  protected:
83  std::vector<double> x_values_;
84  std::vector<T> y_values_;
85  mutable std::vector<T> x_values_reversed_;
86  mutable std::vector<T> y_values_reversed_;
87  };
88 
89 
90  // A utility function
97  template <typename FI>
98  bool isNondecreasing(const FI beg, const FI end)
99  {
100  if (beg == end) return true;
101  FI it = beg;
102  ++it;
103  FI prev = beg;
104  for (; it != end; ++it, ++prev) {
105  if (*it < *prev) {
106  return false;
107  }
108  }
109  return true;
110  }
111 
112 
113 
114  // Member implementations.
115 
116  template<typename T>
117  inline
120  {
121  }
122 
123  template<typename T>
124  inline
126  ::NonuniformTableLinear(const std::vector<double>& x_values,
127  const std::vector<T>& y_values)
128  : x_values_(x_values), y_values_(y_values)
129  {
130  assert(isNondecreasing(x_values.begin(), x_values.end()));
131  }
132 
133  template<typename T>
134  inline std::pair<double, double>
137  {
138  return std::make_pair(x_values_[0], x_values_.back());
139  }
140 
141  template<typename T>
142  inline void
144  ::rescaleDomain(std::pair<double, double> new_domain)
145  {
146  const double a = x_values_[0];
147  const double b = x_values_.back();
148  const double c = new_domain.first;
149  const double d = new_domain.second;
150  // x in [a, b] -> x in [c, d]
151  for (int i = 0; i < int(x_values_.size()); ++i) {
152  x_values_[i] = (x_values_[i] - a)*(d - c)/(b - a) + c;
153  }
154  }
155 
156  template<typename T>
157  inline double
159  ::operator()(const double x) const
160  {
161  return Opm::linearInterpolation(x_values_, y_values_, x);
162  }
163 
164  template<typename T>
165  inline double
167  ::derivative(const double x) const
168  {
169  return Opm::linearInterpolationDerivative(x_values_, y_values_, x);
170  }
171 
172  template<typename T>
173  inline double
175  ::inverse(const double y) const
176  {
177  if (y_values_.front() < y_values_.back()) {
178  return Opm::linearInterpolation(y_values_, x_values_, y);
179  } else {
180  if (y_values_reversed_.empty()) {
181  y_values_reversed_ = y_values_;
182  std::reverse(y_values_reversed_.begin(), y_values_reversed_.end());
183  assert(isNondecreasing(y_values_reversed_.begin(), y_values_reversed_.end()));
184  x_values_reversed_ = x_values_;
185  std::reverse(x_values_reversed_.begin(), x_values_reversed_.end());
186  }
187  return Opm::linearInterpolation(y_values_reversed_, x_values_reversed_, y);
188  }
189  }
190 
191  template<typename T>
192  inline bool
195  {
196  return x_values_ == other.x_values_
197  && y_values_ == other.y_values_;
198  }
199 
200 } // namespace Opm
201 
202 #endif // OPM_NONUNIFORMTABLELINEAR_HEADER_INCLUDED
double operator()(const double x) const
Evaluate the value at x.
Definition: NonuniformTableLinear.hpp:159
double linearInterpolation(const std::vector< double > &xv, const std::vector< double > &yv, double x)
Definition: linearInterpolation.hpp:64
double inverse(const double y) const
Evaluate the inverse at y. Requires T to be a double.
Definition: NonuniformTableLinear.hpp:175
This class uses linear interpolation to compute the value (and its derivative) of a function f sample...
Definition: NonuniformTableLinear.hpp:42
Definition: AnisotropicEikonal.hpp:43
bool isNondecreasing(const FI beg, const FI end)
Detect if a sequence is nondecreasing.
Definition: NonuniformTableLinear.hpp:98
std::pair< double, double > domain()
Get the domain.
Definition: NonuniformTableLinear.hpp:136
bool operator==(const NonuniformTableLinear &other) const
Equality operator.
Definition: NonuniformTableLinear.hpp:194
std::vector< T > y_values_reversed_
Definition: NonuniformTableLinear.hpp:86
std::vector< double > x_values_
Definition: NonuniformTableLinear.hpp:83
std::vector< T > y_values_
Definition: NonuniformTableLinear.hpp:84
std::vector< T > x_values_reversed_
Definition: NonuniformTableLinear.hpp:85
NonuniformTableLinear()
Default constructor.
Definition: NonuniformTableLinear.hpp:119
void rescaleDomain(std::pair< double, double > new_domain)
Rescale the domain.
Definition: NonuniformTableLinear.hpp:144
double derivative(const double x) const
Evaluate the derivative at x.
Definition: NonuniformTableLinear.hpp:167
double linearInterpolationDerivative(const std::vector< double > &xv, const std::vector< double > &yv, double x)
Definition: linearInterpolation.hpp:55