28#ifndef OPM_UNIFORM_TABULATED_2D_FUNCTION_HPP
29#define OPM_UNIFORM_TABULATED_2D_FUNCTION_HPP
31#include <opm/common/OpmLog/OpmLog.hpp>
35#include <opm/common/utility/gpuDecorators.hpp>
37#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
45template <
class Scalar,
template<
class>
class Storage = Opm::VectorWithDefaultAllocator>
52 template <
class ScalarT>
66template <
class Scalar,
template<
class>
class Storage = Opm::VectorWithDefaultAllocator>
67class UniformTabulated2DFunction
69 using ContainerT = Storage<Scalar>;
71 OPM_HOST_DEVICE UniformTabulated2DFunction()
76 UniformTabulated2DFunction(Scalar minX, Scalar maxX,
unsigned m,
77 Scalar minY, Scalar maxY,
unsigned n,
79 : samples_(std::move(
samples)), m_(m), n_(n), xMin_(minX), yMin_(minY), xMax_(maxX), yMax_(maxY){
87 Scalar minY, Scalar maxY,
unsigned n)
89 resize(minX, maxX, m, minY, maxY, n);
93 Scalar minY, Scalar maxY,
unsigned n,
94 const std::vector<std::vector<Scalar>>& vals)
96 resize(minX, maxX, m, minY, maxY, n);
98 for (
unsigned i = 0; i < m; ++i)
99 for (
unsigned j = 0; j < n; ++j)
106 Scalar minY, Scalar maxY,
unsigned n,
107 const double vals[200][500])
109 resize(minX, maxX, m, minY, maxY, n);
111 for (
unsigned i = 0; i < m; ++i)
112 for (
unsigned j = 0; j < n; ++j)
119 void resize(Scalar minX, Scalar maxX,
unsigned m,
120 Scalar minY, Scalar maxY,
unsigned n)
122 samples_.resize(m*n);
137 OPM_HOST_DEVICE Scalar
xMin()
const
143 OPM_HOST_DEVICE Scalar
xMax()
const
149 OPM_HOST_DEVICE Scalar
yMin()
const
155 OPM_HOST_DEVICE Scalar
yMax()
const
161 OPM_HOST_DEVICE
unsigned numX()
const
167 OPM_HOST_DEVICE
unsigned numY()
const
173 OPM_HOST_DEVICE
const ContainerT&
samples()
const
180 OPM_HOST_DEVICE Scalar
iToX(
unsigned i)
const
190 OPM_HOST_DEVICE Scalar
jToY(
unsigned j)
const
205 template <
class Evaluation>
206 OPM_HOST_DEVICE Evaluation
xToI(
const Evaluation& x)
const
217 template <
class Evaluation>
218 OPM_HOST_DEVICE Evaluation
yToJ(
const Evaluation& y)
const
224 template <
class Evaluation>
225 OPM_HOST_DEVICE
bool applies(
const Evaluation& x,
const Evaluation& y)
const
240 template <
class Evaluation>
241 OPM_HOST_DEVICE Evaluation
eval(
const Evaluation& x,
243 [[maybe_unused]]
bool extrapolate)
const
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());
259 Evaluation alpha =
xToI(x);
260 Evaluation beta =
yToJ(y);
263 static_cast<unsigned>(
264 std::max(0, std::min(
static_cast<int>(
numX()) - 2,
265 static_cast<int>(scalarValue(alpha)))));
267 static_cast<unsigned>(
268 std::max(0, std::min(
static_cast<int>(
numY()) - 2,
269 static_cast<int>(scalarValue(beta)))));
277 return s1*(1.0 - beta) + s2*beta;
290 return samples_[j*m_ + i];
303 samples_[j*m_ + i] = value;
308 return samples_ == data.samples_ &&
311 xMin_ == data.xMin_ &&
312 xMax_ == data.xMax_ &&
313 yMin_ == data.yMin_ &&
319 template <
class ScalarT>
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()));
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_));
Provides the OPM specific exception classes.
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, GPUContainerType > copy_to_gpu(const PiecewiseLinearTwoPhaseMaterialParams< TraitsT > ¶ms)
Move a PiecewiseLinearTwoPhaseMaterialParams-object to the GPU.
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:285
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ViewType > make_view(PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ContainerType > ¶ms)
this function is intented to make a GPU friendly view of the PiecewiseLinearTwoPhaseMaterialParams
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:312
Definition Exceptions.hpp:40
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