opm-common
Loading...
Searching...
No Matches
UniformTabulated2DFunction.hpp
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*/
28#ifndef OPM_UNIFORM_TABULATED_2D_FUNCTION_HPP
29#define OPM_UNIFORM_TABULATED_2D_FUNCTION_HPP
30
31#include <opm/common/OpmLog/OpmLog.hpp>
33
35#include <opm/common/utility/gpuDecorators.hpp>
37#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
38
39#include <algorithm>
40#include <cassert>
41#include <string>
42#include <vector>
43
44// forward declaration of the class so the function in the next namespace can be declared
45template <class Scalar, template<class> class Storage = Opm::VectorWithDefaultAllocator>
47
48#if HAVE_CUDA
49// declaration of make_view in correct namespace so friend function can be declared in the class
50namespace Opm::gpuistl
51{
52 template <class ScalarT>
54} // namespace Opm::gpuistl
55#endif // HAVE_CUDA
56
57namespace Opm {
58
66template <class Scalar, template<class> class Storage = Opm::VectorWithDefaultAllocator>
67class UniformTabulated2DFunction
68{
69 using ContainerT = Storage<Scalar>;
70public:
71 OPM_HOST_DEVICE UniformTabulated2DFunction()
72 { }
73
74
75 // Intended for construction of UniformTabulated2DFunction<double, GPUBuffer>
76 UniformTabulated2DFunction(Scalar minX, Scalar maxX, unsigned m,
77 Scalar minY, Scalar maxY, unsigned n,
78 ContainerT&& samples)
79 : samples_(std::move(samples)), m_(m), n_(n), xMin_(minX), yMin_(minY), xMax_(maxX), yMax_(maxY){
80 }
81
86 UniformTabulated2DFunction(Scalar minX, Scalar maxX, unsigned m,
87 Scalar minY, Scalar maxY, unsigned n)
88 {
89 resize(minX, maxX, m, minY, maxY, n);
90 }
91
92 UniformTabulated2DFunction(Scalar minX, Scalar maxX, unsigned m,
93 Scalar minY, Scalar maxY, unsigned n,
94 const std::vector<std::vector<Scalar>>& vals)
95 {
96 resize(minX, maxX, m, minY, maxY, n);
97
98 for (unsigned i = 0; i < m; ++i)
99 for (unsigned j = 0; j < n; ++j)
100 this->setSamplePoint(i, j, vals[i][j]);
101 }
102
103 // Both CO2Tables and H2Tables have values of dimes [200][500]
104 // suboptimal hardcoding for now but easier than templating this function etc
105 UniformTabulated2DFunction(Scalar minX, Scalar maxX, unsigned m,
106 Scalar minY, Scalar maxY, unsigned n,
107 const double vals[200][500])
108 {
109 resize(minX, maxX, m, minY, maxY, n);
110
111 for (unsigned i = 0; i < m; ++i)
112 for (unsigned j = 0; j < n; ++j)
113 this->setSamplePoint(i, j, vals[i][j]);
114 }
115
119 void resize(Scalar minX, Scalar maxX, unsigned m,
120 Scalar minY, Scalar maxY, unsigned n)
121 {
122 samples_.resize(m*n);
123
124 m_ = m;
125 n_ = n;
126
127 xMin_ = minX;
128 xMax_ = maxX;
129
130 yMin_ = minY;
131 yMax_ = maxY;
132 }
133
137 OPM_HOST_DEVICE Scalar xMin() const
138 { return xMin_; }
139
143 OPM_HOST_DEVICE Scalar xMax() const
144 { return xMax_; }
145
149 OPM_HOST_DEVICE Scalar yMin() const
150 { return yMin_; }
151
155 OPM_HOST_DEVICE Scalar yMax() const
156 { return yMax_; }
157
161 OPM_HOST_DEVICE unsigned numX() const
162 { return m_; }
163
167 OPM_HOST_DEVICE unsigned numY() const
168 { return n_; }
169
173 OPM_HOST_DEVICE const ContainerT& samples() const
174 { return samples_; }
175
176
180 OPM_HOST_DEVICE Scalar iToX(unsigned i) const
181 {
182 assert(i < numX());
183
184 return xMin() + i*(xMax() - xMin())/(numX() - 1);
185 }
186
190 OPM_HOST_DEVICE Scalar jToY(unsigned j) const
191 {
192 assert(j < numY());
193
194 return yMin() + j*(yMax() - yMin())/(numY() - 1);
195 }
196
205 template <class Evaluation>
206 OPM_HOST_DEVICE Evaluation xToI(const Evaluation& x) const
207 { return (x - xMin())/(xMax() - xMin())*(numX() - 1); }
208
217 template <class Evaluation>
218 OPM_HOST_DEVICE Evaluation yToJ(const Evaluation& y) const
219 { return (y - yMin())/(yMax() - yMin())*(numY() - 1); }
220
224 template <class Evaluation>
225 OPM_HOST_DEVICE bool applies(const Evaluation& x, const Evaluation& y) const
226 {
227 return
228 xMin() <= x && x <= xMax() &&
229 yMin() <= y && y <= yMax();
230 }
231
240 template <class Evaluation>
241 OPM_HOST_DEVICE Evaluation eval(const Evaluation& x,
242 const Evaluation& y,
243 [[maybe_unused]] bool extrapolate) const
244 {
245#ifndef NDEBUG
246#if !OPM_IS_INSIDE_DEVICE_FUNCTION
247 if (!extrapolate && !applies(x,y)) {
248 std::string msg = "Attempt to get tabulated value for ("
249 +std::to_string(double(scalarValue(x)))+", "+std::to_string(double(scalarValue(y)))
250 +") on a table of extent "
251 +std::to_string(xMin())+" to "+std::to_string(xMax())+" times "
252 +std::to_string(yMin())+" to "+std::to_string(yMax());
253
254 throw NumericalProblem(msg);
255 }
256#endif
257#endif
258
259 Evaluation alpha = xToI(x);
260 Evaluation beta = yToJ(y);
261
262 unsigned i =
263 static_cast<unsigned>(
264 std::max(0, std::min(static_cast<int>(numX()) - 2,
265 static_cast<int>(scalarValue(alpha)))));
266 unsigned j =
267 static_cast<unsigned>(
268 std::max(0, std::min(static_cast<int>(numY()) - 2,
269 static_cast<int>(scalarValue(beta)))));
270
271 alpha -= i;
272 beta -= j;
273
274 // bi-linear interpolation
275 const Evaluation& s1 = getSamplePoint(i, j)*(1.0 - alpha) + getSamplePoint(i + 1, j)*alpha;
276 const Evaluation& s2 = getSamplePoint(i, j + 1)*(1.0 - alpha) + getSamplePoint(i + 1, j + 1)*alpha;
277 return s1*(1.0 - beta) + s2*beta;
278 }
279
285 OPM_HOST_DEVICE Scalar getSamplePoint(unsigned i, unsigned j) const
286 {
287 assert(i < m_);
288 assert(j < n_);
289
290 return samples_[j*m_ + i];
291 }
292
298 void setSamplePoint(unsigned i, unsigned j, Scalar value)
299 {
300 assert(i < m_);
301 assert(j < n_);
302
303 samples_[j*m_ + i] = value;
304 }
305
306 OPM_HOST_DEVICE bool operator==(const UniformTabulated2DFunction<Scalar>& data) const
307 {
308 return samples_ == data.samples_ &&
309 m_ == data.m_ &&
310 n_ == data.n_ &&
311 xMin_ == data.xMin_ &&
312 xMax_ == data.xMax_ &&
313 yMin_ == data.yMin_ &&
314 yMax_ == data.yMax_;
315 }
316
317private:
318 #if HAVE_CUDA
319 template <class ScalarT>
321 #endif
322
323 // the vector which contains the values of the sample points
324 // f(x_i, y_j). don't use this directly, use getSamplePoint(i,j)
325 // instead!
326 ContainerT samples_;
327
328 // the number of sample points in x direction
329 unsigned m_;
330
331 // the number of sample points in y direction
332 unsigned n_;
333
334 // the range of the tabulation on the x axis
335 Scalar xMin_;
336 Scalar xMax_;
337
338 // the range of the tabulation on the y axis
339 Scalar yMin_;
340 Scalar yMax_;
341};
342
343} // namespace Opm
344
345#if HAVE_CUDA
346namespace Opm::gpuistl {
347 template<class ScalarT>
348 UniformTabulated2DFunction<ScalarT, GpuBuffer>
349 copy_to_gpu(const UniformTabulated2DFunction<ScalarT>& tab) {
350 return UniformTabulated2DFunction<ScalarT, GpuBuffer>(tab.xMin(), tab.xMax(), tab.numX(), tab.yMin(), tab.yMax(), tab.numY(), GpuBuffer(tab.samples()));
351 }
352
353 template <class ScalarT>
354 UniformTabulated2DFunction<ScalarT, GpuView>
355 make_view(UniformTabulated2DFunction<ScalarT, GpuBuffer>& tab) {
356 return UniformTabulated2DFunction<ScalarT, GpuView>(tab.xMin(), tab.xMax(), tab.numX(), tab.yMin(), tab.yMax(), tab.numY(), make_view(tab.samples_));
357 }
358} // namespace Opm::gpuistl
359
360#endif // HAVE_CUDA
361
362#endif
Provides the OPM specific exception classes.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, GPUContainerType > copy_to_gpu(const PiecewiseLinearTwoPhaseMaterialParams< TraitsT > &params)
Move a PiecewiseLinearTwoPhaseMaterialParams-object to the GPU.
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:285
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ViewType > make_view(PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ContainerType > &params)
this function is intented to make a GPU friendly view of the PiecewiseLinearTwoPhaseMaterialParams
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:312
Definition Exceptions.hpp:40
Implements a scalar function that depends on two variables and which is sampled on an uniform X-Y gri...
Definition UniformTabulated2DFunction.hpp:68
void resize(Scalar minX, Scalar maxX, unsigned m, Scalar minY, Scalar maxY, unsigned n)
Resize the tabulation to a new range.
Definition UniformTabulated2DFunction.hpp:119
OPM_HOST_DEVICE bool applies(const Evaluation &x, const Evaluation &y) const
Returns true iff a coordinate lies in the tabulated range.
Definition UniformTabulated2DFunction.hpp:225
OPM_HOST_DEVICE Scalar iToX(unsigned i) const
Return the position on the x-axis of the i-th interval.
Definition UniformTabulated2DFunction.hpp:180
OPM_HOST_DEVICE unsigned numY() const
Returns the number of sampling points in Y direction.
Definition UniformTabulated2DFunction.hpp:167
OPM_HOST_DEVICE Scalar yMin() const
Returns the minimum of the Y coordinate of the sampling points.
Definition UniformTabulated2DFunction.hpp:149
OPM_HOST_DEVICE Scalar xMin() const
Returns the minimum of the X coordinate of the sampling points.
Definition UniformTabulated2DFunction.hpp:137
OPM_HOST_DEVICE Evaluation yToJ(const Evaluation &y) const
Return the interval index of a given position on the y-axis.
Definition UniformTabulated2DFunction.hpp:218
OPM_HOST_DEVICE Scalar getSamplePoint(unsigned i, unsigned j) const
Get the value of the sample point which is at the intersection of the -th interval of the x-Axis and ...
Definition UniformTabulated2DFunction.hpp:285
OPM_HOST_DEVICE const ContainerT & samples() const
Returns the sampling points.
Definition UniformTabulated2DFunction.hpp:173
void setSamplePoint(unsigned i, unsigned j, Scalar value)
Set the value of the sample point which is at the intersection of the -th interval of the x-Axis and ...
Definition UniformTabulated2DFunction.hpp:298
OPM_HOST_DEVICE Scalar yMax() const
Returns the maximum of the Y coordinate of the sampling points.
Definition UniformTabulated2DFunction.hpp:155
UniformTabulated2DFunction(Scalar minX, Scalar maxX, unsigned m, Scalar minY, Scalar maxY, unsigned n)
Constructor where the tabulation parameters are already provided.
Definition UniformTabulated2DFunction.hpp:86
OPM_HOST_DEVICE Scalar xMax() const
Returns the maximum of the X coordinate of the sampling points.
Definition UniformTabulated2DFunction.hpp:143
OPM_HOST_DEVICE Evaluation eval(const Evaluation &x, const Evaluation &y, bool extrapolate) const
Evaluate the function at a given (x,y) position.
Definition UniformTabulated2DFunction.hpp:241
OPM_HOST_DEVICE Scalar jToY(unsigned j) const
Return the position on the y-axis of the j-th interval.
Definition UniformTabulated2DFunction.hpp:190
OPM_HOST_DEVICE unsigned numX() const
Returns the number of sampling points in X direction.
Definition UniformTabulated2DFunction.hpp:161
OPM_HOST_DEVICE Evaluation xToI(const Evaluation &x) const
Return the interval index of a given position on the x-axis.
Definition UniformTabulated2DFunction.hpp:206
Definition UniformTabulated2DFunction.hpp:46
Convience header to include the gpuistl headers if HAVE_CUDA is defined.
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30