version 3.10-dev
spline.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_SPLINE_HH
13#define DUMUX_SPLINE_HH
14
15#include "fixedlengthspline_.hh"
16#include "variablelengthspline_.hh"
17#include "splinecommon_.hh"
18
19namespace Dumux {
20
41template<class Scalar, int numSamples = 2>
42class Spline : public FixedLengthSpline_<Scalar, numSamples>
43{
44public:
51 { };
52
59 template <class ScalarArray>
60 Spline(const ScalarArray &x,
61 const ScalarArray &y)
62 { this->setXYArrays(numSamples, x, y); }
63
69 template <class PointArray>
70 Spline(const PointArray &points)
71 { this->setArrayOfPoints(numSamples, points); }
72
81 template <class ScalarArray>
82 Spline(const ScalarArray &x,
83 const ScalarArray &y,
84 Scalar m0,
85 Scalar m1)
86 { this->setXYArrays(numSamples, x, y, m0, m1); }
87
95 template <class PointArray>
96 Spline(const PointArray &points,
97 Scalar m0,
98 Scalar m1)
99 { this->setArrayOfPoints(numSamples, points, m0, m1); }
100};
101
123template<class Scalar>
124class Spline<Scalar, /*numSamples=*/-1> : public VariableLengthSpline_<Scalar>
125{
126public:
133 { }
134
142 template <class ScalarArrayX, class ScalarArrayY>
143 Spline(int nSamples,
144 const ScalarArrayX &x,
145 const ScalarArrayY &y)
146 { this->setXYArrays(nSamples, x, y); }
147
154 template <class PointArray>
155 Spline(int nSamples,
156 const PointArray &points)
157 { this->setArrayOfPoints(nSamples, points); }
158
165 template <class ScalarContainer>
166 Spline(const ScalarContainer &x,
167 const ScalarContainer &y)
168 { this->setXYContainers(x, y); }
169
175 template <class PointContainer>
176 Spline(const PointContainer &points)
177 { this->setContainerOfPoints(points); }
178
188 template <class ScalarArray>
189 Spline(int nSamples,
190 const ScalarArray &x,
191 const ScalarArray &y,
192 Scalar m0,
193 Scalar m1)
194 { this->setXYArrays(nSamples, x, y, m0, m1); }
195
204 template <class PointArray>
205 Spline(int nSamples,
206 const PointArray &points,
207 Scalar m0,
208 Scalar m1)
209 { this->setArrayOfPoints(nSamples, points, m0, m1); }
210
219 template <class ScalarContainerX, class ScalarContainerY>
220 Spline(const ScalarContainerX &x,
221 const ScalarContainerY &y,
222 Scalar m0,
223 Scalar m1)
224 { this->setXYContainers(x, y, m0, m1); }
225
233 template <class PointContainer>
234 Spline(const PointContainer &points,
235 Scalar m0,
236 Scalar m1)
237 { this->setContainerOfPoints(points, m0, m1); }
238};
239
243template<class Scalar>
244class Spline<Scalar, /*numSamples=*/0>
245// Splines with zero sampling points do not make sense!
246{ private: Spline() { }; };
247
251template<class Scalar>
252class Spline<Scalar, /*numSamples=*/1>
253// Splines with one sampling point do not make sense!
254{ private: Spline() { }; };
255
262template<class Scalar>
263class Spline<Scalar, 2> : public SplineCommon_<Scalar, Spline<Scalar, 2> >
264{
265 friend class SplineCommon_<Scalar, Spline<Scalar, 2> >;
266 using Vector = Dune::FieldVector<Scalar, 2>;
267 using Matrix = Dune::FieldMatrix<Scalar, 2, 2>;
268
269public:
271 {};
272
281 template <class ScalarArrayX, class ScalarArrayY>
282 Spline(const ScalarArrayX &x,
283 const ScalarArrayY &y,
284 Scalar m0, Scalar m1)
285 { setXYArrays(2, x, y, m0, m1); }
286
294 template <class PointArray>
295 Spline(const PointArray &points,
296 Scalar m0,
297 Scalar m1)
298 { this->setArrayOfPoints(2, points, m0, m1); }
299
310 Spline(Scalar x0, Scalar x1,
311 Scalar y0, Scalar y1,
312 Scalar m0, Scalar m1)
313 {
314 set(x0, x1,
315 y0, y1,
316 m0, m1);
317 };
318
322 int numSamples() const
323 { return 2; }
324
336 void set(Scalar x0, Scalar x1,
337 Scalar y0, Scalar y1,
338 Scalar m0, Scalar m1)
339 {
340 Matrix M(numSamples());
341 Vector d;
342 assignXY_(x0, x1, y0, y1);
343 this->makeFullSystem_(M, d, m0, m1);
344
345 // solve for the moments
346 M.solve(m_, d);
347 }
348
359 template <class ScalarContainer>
360 void setXYArrays(int nSamples,
361 const ScalarContainer &x,
362 const ScalarContainer &y,
363 Scalar m0, Scalar m1)
364 {
365 assert(nSamples == 2);
366 set(x[0], x[1], y[0], y[1], m0, m1);
367 }
368
378 template <class ScalarContainerX, class ScalarContainerY>
379 void setXYContainers(const ScalarContainerX &x,
380 const ScalarContainerY &y,
381 Scalar m0, Scalar m1)
382 {
383 assert(x.size() == y.size());
384 assert(x.size() == 2);
385
386 Matrix M(numSamples());
387 Vector d;
388
389 typename ScalarContainerX::const_iterator xIt0 = x.begin();
390 typename ScalarContainerX::const_iterator xIt1 = xIt0;
391 ++xIt1;
392 typename ScalarContainerY::const_iterator yIt0 = y.begin();
393 typename ScalarContainerY::const_iterator yIt1 = yIt0;
394 ++yIt1;
395 set(*xIt0, *xIt1, *yIt0, *yIt1);
396 }
397
407 template <class PointArray>
408 void setArrayOfPoints(int nSamples,
409 const PointArray &points,
410 Scalar m0,
411 Scalar m1)
412 {
413 assert(nSamples == 2);
414
415 set(points[0][0],
416 points[1][0],
417 points[0][1],
418 points[1][1],
419 m0, m1);
420 }
421
430 template <class PointContainer>
431 void setContainerOfPoints(const PointContainer &points,
432 Scalar m0,
433 Scalar m1)
434 {
435 assert(points.size() == 2);
436
437 Matrix M;
438 Vector d;
439 typename PointContainer::const_iterator it0 = points.begin();
440 typename PointContainer::const_iterator it1 = it0;
441 ++it1;
442
443 set((*it0)[0],
444 (*it0)[1],
445 (*it1)[0],
446 (*it1)[1],
447 m0, m1);
448 }
449
458 template <class TupleContainer>
459 void setContainerOfTuples(const TupleContainer &tuples,
460 Scalar m0,
461 Scalar m1)
462 {
463 assert(tuples.size() == 2);
464
465 typename TupleContainer::const_iterator it0 = tuples.begin();
466 typename TupleContainer::const_iterator it1 = it0;
467 ++it1;
468
469 set(std::get<0>(*it0),
470 std::get<1>(*it0),
471 std::get<0>(*it1),
472 std::get<1>(*it1),
473 m0, m1);
474 }
475
476protected:
477 void assignXY_(Scalar x0, Scalar x1,
478 Scalar y0, Scalar y1)
479 {
480 if (x0 > x1) {
481 xPos_[0] = x1;
482 xPos_[1] = x0;
483 yPos_[0] = y1;
484 yPos_[1] = y0;
485 }
486 else {
487 xPos_[0] = x0;
488 xPos_[1] = x1;
489 yPos_[0] = y0;
490 yPos_[1] = y1;
491 }
492 };
493
497 Scalar x_(int i) const
498 { return xPos_[i]; }
499
503 Scalar y_(int i) const
504 { return yPos_[i]; }
505
510 Scalar moment_(int i) const
511 { return m_[i]; }
512
513 Vector xPos_;
514 Vector yPos_;
515 Vector m_;
516};
517
518} // end namespace Dumux
519
520#endif
Spline(const ScalarContainer &x, const ScalarContainer &y)
Convenience constructor for a natural spline.
Definition: spline.hh:166
Spline()
Default constructor for a spline.
Definition: spline.hh:132
Spline(int nSamples, const ScalarArray &x, const ScalarArray &y, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:189
Spline(const PointContainer &points, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:234
Spline(int nSamples, const PointArray &points)
Convenience constructor for a natural spline.
Definition: spline.hh:155
Spline(const PointContainer &points)
Convenience constructor for a natural spline.
Definition: spline.hh:176
Spline(int nSamples, const PointArray &points, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:205
Spline(int nSamples, const ScalarArrayX &x, const ScalarArrayY &y)
Convenience constructor for a natural spline.
Definition: spline.hh:143
Spline(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:220
int numSamples() const
Returns the number of sampling points.
Definition: spline.hh:322
void set(Scalar x0, Scalar x1, Scalar y0, Scalar y1, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of the spline.
Definition: spline.hh:336
Scalar y_(int i) const
Returns the y coordinate of the i-th sampling point.
Definition: spline.hh:503
void setContainerOfPoints(const PointContainer &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes from an STL-like container of points.
Definition: spline.hh:431
Vector m_
Definition: spline.hh:515
Spline(Scalar x0, Scalar x1, Scalar y0, Scalar y1, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:310
Spline(const PointArray &points, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:295
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of the spline.
Definition: spline.hh:379
Vector yPos_
Definition: spline.hh:514
void setContainerOfTuples(const TupleContainer &tuples, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes from an STL-like container of tuples.
Definition: spline.hh:459
Spline(const ScalarArrayX &x, const ScalarArrayY &y, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:282
void assignXY_(Scalar x0, Scalar x1, Scalar y0, Scalar y1)
Definition: spline.hh:477
Scalar moment_(int i) const
Returns the moment (i.e. second derivative) of the spline at the i-th sampling point.
Definition: spline.hh:510
Spline()
Definition: spline.hh:270
void setXYArrays(int nSamples, const ScalarContainer &x, const ScalarContainer &y, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of the spline.
Definition: spline.hh:360
Scalar x_(int i) const
Returns the x coordinate of the i-th sampling point.
Definition: spline.hh:497
void setArrayOfPoints(int nSamples, const PointArray &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of the spline.
Definition: spline.hh:408
Vector xPos_
Definition: spline.hh:513
A 3rd order polynomial spline.
Definition: spline.hh:43
Spline()
Default constructor for a spline.
Definition: spline.hh:50
Spline(const PointArray &points, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:96
Spline(const PointArray &points)
Convenience constructor for a full spline.
Definition: spline.hh:70
Spline(const ScalarArray &x, const ScalarArray &y)
Convenience constructor for a full spline.
Definition: spline.hh:60
Spline(const ScalarArray &x, const ScalarArray &y, Scalar m0, Scalar m1)
Convenience constructor for a full spline.
Definition: spline.hh:82
Definition: adapt.hh:17