opm-common
Loading...
Searching...
No Matches
Tabulated1DFunction.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*/
27#ifndef OPM_TABULATED_1D_FUNCTION_HPP
28#define OPM_TABULATED_1D_FUNCTION_HPP
29
30#include <opm/common/OpmLog/OpmLog.hpp>
32
33#include <algorithm>
34#include <cassert>
35#include <iosfwd>
36#include <stdexcept>
37#include <vector>
38
39namespace Opm {
40
42 size_t value;
43};
44
49template <class Scalar>
51{
52public:
60
69 template <class ScalarArrayX, class ScalarArrayY>
70 Tabulated1DFunction(size_t nSamples,
71 const ScalarArrayX& x,
72 const ScalarArrayY& y,
73 bool sortInputs = true)
74 { this->setXYArrays(nSamples, x, y, sortInputs); }
75
85 template <class ScalarContainer>
86 Tabulated1DFunction(const ScalarContainer& x,
87 const ScalarContainer& y,
88 bool sortInputs = true)
89 { this->setXYContainers(x, y, sortInputs); }
90
98 template <class PointContainer>
99 explicit Tabulated1DFunction(const PointContainer& points,
100 bool sortInputs = true)
101 { this->setContainerOfTuples(points, sortInputs); }
102
108 template <class ScalarArrayX, class ScalarArrayY>
109 void setXYArrays(size_t nSamples,
110 const ScalarArrayX& x,
111 const ScalarArrayY& y,
112 bool sortInputs = true)
113 {
114 assert(nSamples > 1);
115
116 resizeArrays_(nSamples);
117 for (size_t i = 0; i < nSamples; ++i) {
118 xValues_[i] = x[i];
119 yValues_[i] = y[i];
120 }
121
122 if (sortInputs)
123 sortInput_();
124 else if (xValues_[0] > xValues_[numSamples() - 1])
125 reverseSamplingPoints_();
126 }
127
133 template <class ScalarContainerX, class ScalarContainerY>
134 void setXYContainers(const ScalarContainerX& x,
135 const ScalarContainerY& y,
136 bool sortInputs = true)
137 {
138 assert(x.size() == y.size());
139
140 resizeArrays_(x.size());
141 if (x.size() > 0) {
142 std::ranges::copy(x, xValues_.begin());
143 std::ranges::copy(y, yValues_.begin());
144
145 if (sortInputs)
146 sortInput_();
147 else if (xValues_[0] > xValues_[numSamples() - 1])
148 reverseSamplingPoints_();
149 }
150 }
151
155 template <class PointArray>
156 void setArrayOfPoints(size_t nSamples,
157 const PointArray& points,
158 bool sortInputs = true)
159 {
160 // a linear function with less than two sampling points? what an incredible
161 // bad idea!
162 assert(nSamples > 1);
163
164 resizeArrays_(nSamples);
165 for (size_t i = 0; i < nSamples; ++i) {
166 xValues_[i] = points[i][0];
167 yValues_[i] = points[i][1];
168 }
169
170 if (sortInputs)
171 sortInput_();
172 else if (xValues_[0] > xValues_[numSamples() - 1])
173 reverseSamplingPoints_();
174 }
175
190 template <class XYContainer>
191 void setContainerOfTuples(const XYContainer& points,
192 bool sortInputs = true)
193 {
194 // a linear function with less than two sampling points? what an incredible
195 // bad idea!
196 assert(points.size() > 1);
197
198 resizeArrays_(points.size());
199 typename XYContainer::const_iterator it = points.begin();
200 typename XYContainer::const_iterator endIt = points.end();
201 for (unsigned i = 0; it != endIt; ++i, ++it) {
202 xValues_[i] = std::get<0>(*it);
203 yValues_[i] = std::get<1>(*it);
204 }
205
206 if (sortInputs)
207 sortInput_();
208 else if (xValues_[0] > xValues_[numSamples() - 1])
209 reverseSamplingPoints_();
210 }
211
215 size_t numSamples() const
216 { return xValues_.size(); }
217
221 Scalar xMin() const
222 { return xValues_[0]; }
223
227 Scalar xMax() const
228 { return xValues_[numSamples() - 1]; }
229
233 Scalar xAt(size_t i) const
234 { return xValues_[i]; }
235
236 const std::vector<Scalar>& xValues() const
237 { return xValues_; }
238
239 const std::vector<Scalar>& yValues() const
240 { return yValues_; }
241
245 Scalar valueAt(size_t i) const
246 { return yValues_[i]; }
247
251 template <class Evaluation>
252 bool applies(const Evaluation& x) const
253 { return xValues_[0] <= x && x <= xValues_[numSamples() - 1]; }
254
264 template <class Evaluation>
265 Evaluation eval(const Evaluation& x, bool extrapolate = false) const
266 {
267 SegmentIndex segIdx = findSegmentIndex(x, extrapolate);
268 return eval(x, segIdx);
269 }
270
271 template <class Evaluation>
272 Evaluation eval(const Evaluation& x, SegmentIndex segIdxIn) const
273 {
274 size_t segIdx = segIdxIn.value;
275 Scalar x0 = xValues_[segIdx];
276 Scalar x1 = xValues_[segIdx + 1];
277
278 Scalar y0 = yValues_[segIdx];
279 Scalar y1 = yValues_[segIdx + 1];
280
281 return y0 + (y1 - y0)*(x - x0)/(x1 - x0);
282 }
283
295 template <class Evaluation>
296 Evaluation evalDerivative(const Evaluation& x, bool extrapolate = false) const
297 {
298 size_t segIdx = findSegmentIndex(x, extrapolate).value;
299 return evalDerivative_(x, segIdx);
300 }
301
316 template <class Evaluation>
317 Evaluation evalSecondDerivative(const Evaluation&, bool = false) const
318 { return 0.0; }
319
334 template <class Evaluation>
335 Evaluation evalThirdDerivative(const Evaluation&, bool = false) const
336 { return 0.0; }
337
346 int monotonic(Scalar x0, Scalar x1,
347 [[maybe_unused]] bool extrapolate = false) const
348 {
349 assert(x0 != x1);
350
351 // make sure that x0 is smaller than x1
352 if (x0 > x1)
353 std::swap(x0, x1);
354
355 assert(x0 < x1);
356
357 int r = 3;
358 if (x0 < xMin()) {
359 assert(extrapolate);
360
361 x0 = xMin();
362 };
363
364 size_t i = findSegmentIndex(x0, extrapolate).value;
365 if (xValues_[i + 1] >= x1) {
366 // interval is fully contained within a single function
367 // segment
368 updateMonotonicity_(i, r);
369 return r;
370 }
371
372 // the first segment overlaps with the specified interval
373 // partially
374 updateMonotonicity_(i, r);
375 ++ i;
376
377 // make sure that the segments which are completly in the
378 // interval [x0, x1] all exhibit the same monotonicity.
379 size_t iEnd = findSegmentIndex(x1, extrapolate).value;
380 for (; i < iEnd - 1; ++i) {
381 updateMonotonicity_(i, r);
382 if (!r)
383 return 0;
384 }
385
386 // if the user asked for a part of the function which is
387 // extrapolated, we need to check the slope at the function's
388 // endpoint
389 if (x1 > xMax()) {
390 assert(extrapolate);
391
392 Scalar m = evalDerivative_(xMax(), /*segmentIdx=*/numSamples() - 2);
393 if (m < 0)
394 return (r < 0 || r==3) ? -1 : 0;
395 else if (m > 0)
396 return r > 0 ? 1 : 0;
397
398 return r;
399 }
400
401 // check for the last segment
402 updateMonotonicity_(iEnd, r);
403
404 return r;
405 }
406
411 int monotonic() const
412 { return monotonic(xMin(), xMax()); }
413
430 void printCSV(Scalar xi0, Scalar xi1, unsigned k, std::ostream& os) const;
431
432 bool operator==(const Tabulated1DFunction<Scalar>& data) const {
433 return xValues_ == data.xValues_ &&
434 yValues_ == data.yValues_;
435 }
436
437 template <class Evaluation>
438 SegmentIndex findSegmentIndex(const Evaluation& x, bool extrapolate = false) const
439 {
440 if (!isfinite(x)) {
441 throw std::runtime_error("We can not search for extrapolation/interpolation "
442 "segment in an 1D table for non-finite value " +
443 std::to_string(getValue(x)) + " .");
444 }
445
446 if (!extrapolate && !applies(x))
447 throw std::logic_error("Trying to evaluate a tabulated function outside of its range");
448
449 // we need at least two sampling points!
450 if (numSamples() < 2) {
451 throw std::logic_error("We need at least two sampling points to "
452 "do interpolation/extrapolation, "
453 "and the table only contains " +
454 std::to_string(numSamples()) +
455 " sampling points");
456 }
457
458 if (x <= xValues_[1])
459 return SegmentIndex{0};
460 else if (x >= xValues_[xValues_.size() - 2])
461 return SegmentIndex{xValues_.size() - 2};
462 else {
463 // bisection
464 size_t lowerIdx = 1;
465 size_t upperIdx = xValues_.size() - 2;
466 while (lowerIdx + 1 < upperIdx) {
467 size_t pivotIdx = (lowerIdx + upperIdx) / 2;
468 if (x < xValues_[pivotIdx])
469 upperIdx = pivotIdx;
470 else
471 lowerIdx = pivotIdx;
472 }
473
474 if (xValues_[lowerIdx] > x || x > xValues_[lowerIdx + 1]) {
475 std::string msg = "Problematic interpolation/extrapolation "
476 "segment is found for the input value " +
477 std::to_string(Opm::getValue(x)) +
478 "\nthe lower index of the found segment is " +
479 std::to_string(lowerIdx) +
480 ", the size of the table is " +
481 std::to_string(numSamples()) +
482 ",\nand the end values of the found segment are " +
483 std::to_string(xValues_[lowerIdx]) +
484 " and " +
485 std::to_string(xValues_[lowerIdx + 1]) +
486 ", respectively.\n";
487 msg += "Outputting the problematic table for more information "
488 "(with *** marking the found segment):";
489 for (size_t i = 0; i < numSamples(); ++i) {
490 if (i % 10 == 0)
491 msg += "\n";
492 if (i == lowerIdx)
493 msg += " ***";
494 msg += " " + std::to_string(xValues_[i]);
495 if (i == lowerIdx + 1)
496 msg += " ***";
497 }
498 msg += "\n";
499 OpmLog::debug(msg);
500 throw std::runtime_error(msg);
501 }
502 return SegmentIndex{lowerIdx};
503 }
504 }
505
506private:
507 template <class Evaluation>
508 Evaluation evalDerivative_(const Evaluation& x, size_t segIdx) const
509 {
510
511 Scalar x0 = xValues_[segIdx];
512 Scalar x1 = xValues_[segIdx + 1];
513
514 Scalar y0 = yValues_[segIdx];
515 Scalar y1 = yValues_[segIdx + 1];
516
517 Evaluation ret = blank(x);
518 ret = (y1 - y0)/(x1 - x0);
519 return ret;
520 }
521
522 // returns the monotonicity of a segment
523 //
524 // The return value have the following meaning:
525 //
526 // 3: function is constant within interval [x0, x1]
527 // 1: function is monotonously increasing in the specified interval
528 // 0: function is not monotonic in the specified interval
529 // -1: function is monotonously decreasing in the specified interval
530 int updateMonotonicity_(size_t i, int& r) const
531 {
532 if (yValues_[i] < yValues_[i + 1]) {
533 // monotonically increasing?
534 if (r == 3 || r == 1)
535 r = 1;
536 else
537 r = 0;
538 return 1;
539 }
540 else if (yValues_[i] > yValues_[i + 1]) {
541 // monotonically decreasing?
542 if (r == 3 || r == -1)
543 r = -1;
544 else
545 r = 0;
546 return -1;
547 }
548
549 return 3;
550 }
551
555 struct ComparatorX_
556 {
557 explicit ComparatorX_(const std::vector<Scalar>& x)
558 : x_(x)
559 {}
560
561 bool operator ()(size_t idxA, size_t idxB) const
562 { return x_.at(idxA) < x_.at(idxB); }
563
564 const std::vector<Scalar>& x_;
565 };
566
570 void sortInput_()
571 {
572 size_t n = numSamples();
573
574 // create a vector containing 0...n-1
575 std::vector<unsigned> idxVector(n);
576 for (unsigned i = 0; i < n; ++i)
577 idxVector[i] = i;
578
579 // sort the indices according to the x values of the sample
580 // points
581 ComparatorX_ cmp(xValues_);
582 std::ranges::sort(idxVector, cmp);
583
584 // reorder the sample points
585 std::vector<Scalar> tmpX(n), tmpY(n);
586 for (size_t i = 0; i < idxVector.size(); ++ i) {
587 tmpX[i] = xValues_[idxVector[i]];
588 tmpY[i] = yValues_[idxVector[i]];
589 }
590 xValues_ = tmpX;
591 yValues_ = tmpY;
592 }
593
598 void reverseSamplingPoints_()
599 {
600 // reverse the arrays
601 size_t n = numSamples();
602 for (size_t i = 0; i <= (n - 1)/2; ++i) {
603 std::swap(xValues_[i], xValues_[n - i - 1]);
604 std::swap(yValues_[i], yValues_[n - i - 1]);
605 }
606 }
607
611 void resizeArrays_(size_t nSamples)
612 {
613 xValues_.resize(nSamples);
614 yValues_.resize(nSamples);
615 }
616
617 std::vector<Scalar> xValues_;
618 std::vector<Scalar> yValues_;
619};
620
621} // namespace Opm
622
623#endif
A number of commonly used algebraic functions for the localized OPM automatic differentiation (AD) fr...
Evaluation eval(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition Tabulated1DFunction.hpp:265
void printCSV(Scalar xi0, Scalar xi1, unsigned k, std::ostream &os) const
Prints k tuples of the format (x, y, dx/dy, isMonotonic) to stdout.
Definition Tabulated1DFunction.cpp:33
size_t numSamples() const
Returns the number of sampling points.
Definition Tabulated1DFunction.hpp:215
Evaluation evalDerivative(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline's derivative at a given position.
Definition Tabulated1DFunction.hpp:296
void setArrayOfPoints(size_t nSamples, const PointArray &points, bool sortInputs=true)
Set the sampling points for the piecewise linear function.
Definition Tabulated1DFunction.hpp:156
Scalar xMax() const
Return the x value of the rightmost sampling point.
Definition Tabulated1DFunction.hpp:227
Tabulated1DFunction(const ScalarContainer &x, const ScalarContainer &y, bool sortInputs=true)
Convenience constructor for a piecewise linear function.
Definition Tabulated1DFunction.hpp:86
Scalar xMin() const
Return the x value of the leftmost sampling point.
Definition Tabulated1DFunction.hpp:221
Evaluation evalThirdDerivative(const Evaluation &, bool=false) const
Evaluate the function's third derivative at a given position.
Definition Tabulated1DFunction.hpp:335
Evaluation evalSecondDerivative(const Evaluation &, bool=false) const
Evaluate the function's second derivative at a given position.
Definition Tabulated1DFunction.hpp:317
Tabulated1DFunction(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, bool sortInputs=true)
Convenience constructor for a piecewise linear function.
Definition Tabulated1DFunction.hpp:70
void setContainerOfTuples(const XYContainer &points, bool sortInputs=true)
Set the sampling points of the piecewise linear function using a STL-compatible container of tuple-li...
Definition Tabulated1DFunction.hpp:191
Scalar xAt(size_t i) const
Return the x value of the a sample point with a given index.
Definition Tabulated1DFunction.hpp:233
Tabulated1DFunction(const PointContainer &points, bool sortInputs=true)
Convenience constructor for a piecewise linear function.
Definition Tabulated1DFunction.hpp:99
Tabulated1DFunction()
Default constructor for a piecewise linear function.
Definition Tabulated1DFunction.hpp:58
int monotonic() const
Same as monotonic(x0, x1), but with the entire range of the function as interval.
Definition Tabulated1DFunction.hpp:411
void setXYArrays(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, bool sortInputs=true)
Set the sampling points for the piecewise linear function.
Definition Tabulated1DFunction.hpp:109
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, bool sortInputs=true)
Set the sampling points for the piecewise linear function.
Definition Tabulated1DFunction.hpp:134
int monotonic(Scalar x0, Scalar x1, bool extrapolate=false) const
Returns 1 if the function is monotonically increasing, -1 if the function is mononously decreasing an...
Definition Tabulated1DFunction.hpp:346
bool applies(const Evaluation &x) const
Return true iff the given x is in range [x1, xn].
Definition Tabulated1DFunction.hpp:252
Scalar valueAt(size_t i) const
Return the value of the a sample point with a given index.
Definition Tabulated1DFunction.hpp:245
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition Tabulated1DFunction.hpp:41