3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fixedlengthspline_.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FIXED_LENGTH_SPLINE_HH
25#define DUMUX_FIXED_LENGTH_SPLINE_HH
26
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
29#include <dune/istl/btdmatrix.hh>
30
31#include "splinecommon_.hh"
32
33namespace Dumux {
34
40template<class ScalarT, int nSamples>
42 : public SplineCommon_<ScalarT,
43 FixedLengthSpline_<ScalarT, nSamples> >
44{
45 friend class SplineCommon_<ScalarT, FixedLengthSpline_<ScalarT, nSamples> >;
46
47 using Scalar = ScalarT;
48 using Vector = Dune::FieldVector<Scalar, nSamples>;
49 using BlockVector = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
50 using BTDMatrix = Dune::BTDMatrix<Dune::FieldMatrix<Scalar, 1, 1> >;
51
52protected:
54 : m_(nSamples)
55 {};
56
57public:
61 int numSamples() const
62 { return nSamples; }
63
67 // Full splines //
71
83 template <class ScalarArrayX, class ScalarArrayY>
84 void setXYArrays(int numberSamples,
85 const ScalarArrayX &x,
86 const ScalarArrayY &y,
87 Scalar m0, Scalar m1)
88 {
89 assert(numberSamples == numSamples());
90
91 for (int i = 0; i < numberSamples; ++i) {
92 xPos_[i] = x[i];
93 yPos_[i] = y[i];
94 }
95 makeFullSpline_(m0, m1);
96 }
97
109 template <class ScalarContainerX, class ScalarContainerY>
110 void setXYContainers(const ScalarContainerX &x,
111 const ScalarContainerY &y,
112 Scalar m0, Scalar m1)
113 {
114 assert(x.size() == y.size());
115 assert(x.size() == numSamples());
116
117 std::copy(x.begin(), x.end(), xPos_.begin());
118 std::copy(y.begin(), y.end(), yPos_.begin());
119 makeFullSpline_(m0, m1);
120 }
121
134 template <class PointArray>
135 void setArrayOfPoints(int numberSamples,
136 const PointArray &points,
137 Scalar m0,
138 Scalar m1)
139 {
140 // a spline with no or just one sampling points? what an
141 // incredible bad idea!
142 assert(numberSamples == numSamples());
143
144 for (int i = 0; i < numberSamples; ++i) {
145 xPos_[i] = points[i][0];
146 yPos_[i] = points[i][1];
147 }
148 makeFullSpline_(m0, m1);
149 }
150
163 template <class XYContainer>
164 void setContainerOfPoints(const XYContainer &points,
165 Scalar m0,
166 Scalar m1)
167 {
168 assert(points.size() == numSamples());
169
170 typename XYContainer::const_iterator it = points.begin();
171 typename XYContainer::const_iterator endIt = points.end();
172 for (int i = 0; it != endIt; ++i, ++it) {
173 xPos_[i] = (*it)[0];
174 yPos_[i] = (*it)[1];
175 }
176
177 // make a full spline
178 makeFullSpline_(m0, m1);
179 }
180
196 template <class XYContainer>
197 void setContainerOfTuples(const XYContainer &points,
198 Scalar m0,
199 Scalar m1)
200 {
201 assert(points.size() == numSamples());
202
203 typename XYContainer::const_iterator it = points.begin();
204 typename XYContainer::const_iterator endIt = points.end();
205 for (int i = 0; it != endIt; ++i, ++it) {
206 xPos_[i] = std::get<0>(*it);
207 yPos_[i] = std::get<1>(*it);
208 }
209
210 // make a full spline
211 makeFullSpline_(m0, m1);
212 }
213
217 // Natural splines //
221
231 template <class ScalarArrayX, class ScalarArrayY>
232 void setXYArrays(int numberSamples,
233 const ScalarArrayX &x,
234 const ScalarArrayY &y)
235 {
236 assert(numberSamples == numSamples());
237
238 for (int i = 0; i < numberSamples; ++i) {
239 xPos_[i] = x[i];
240 yPos_[i] = y[i];
241 }
242
244 }
245
257 template <class ScalarContainerX, class ScalarContainerY>
258 void setXYContainers(const ScalarContainerX &x,
259 const ScalarContainerY &y)
260 {
261 assert(x.size() == y.size());
262 assert(x.size() > 1);
263
264 setNumSamples_(x.size());
265 std::copy(x.begin(), x.end(), xPos_.begin());
266 std::copy(y.begin(), y.end(), yPos_.begin());
268 }
269
270
283 template <class PointArray>
284 void setArrayOfPoints(int numberSamples,
285 const PointArray &points)
286 {
287 assert(numberSamples == numSamples());
288
289 for (int i = 0; i < numberSamples; ++i) {
290 xPos_[i] = points[i][0];
291 yPos_[i] = points[i][1];
292 }
294 }
295
307 template <class XYContainer>
308 void setContainerOfPoints(const XYContainer &points)
309 {
310 assert(points.size() == numSamples());
311
312 typename XYContainer::const_iterator it = points.begin();
313 typename XYContainer::const_iterator endIt = points.end();
314 for (int i = 0; it != endIt; ++ i, ++it) {
315 xPos_[i] = (*it)[0];
316 yPos_[i] = (*it)[1];
317 }
318
319 // make a natural spline
321 }
322
337 template <class XYContainer>
338 void setContainerOfTuples(const XYContainer &points)
339 {
340 assert(points.size() == numSamples());
341
342 typename XYContainer::const_iterator it = points.begin();
343 typename XYContainer::const_iterator endIt = points.end();
344 for (int i = 0; it != endIt; ++i, ++it) {
345 xPos_[i] = std::get<0>(*it);
346 yPos_[i] = std::get<1>(*it);
347 }
348
349 // make a natural spline
351 }
352
353protected:
359 void makeFullSpline_(Scalar m0, Scalar m1)
360 {
361 BTDMatrix M(numSamples());
362 BlockVector d(numSamples());
363
364 // create linear system of equations
365 this->makeFullSystem_(M, d, m0, m1);
366
367 // solve for the moments
368 M.solve(m_, d);
369 }
370
377 {
378 BTDMatrix M(numSamples());
379 BlockVector d(numSamples());
380
381 // create linear system of equations
382 this->makeNaturalSystem_(M, d);
383
384 // solve for the moments
385 M.solve(m_, d);
386 }
387
391 Scalar x_(int i) const
392 { return xPos_[i]; }
393
397 Scalar y_(int i) const
398 { return yPos_[i]; }
399
404 Scalar moment_(int i) const
405 { return m_[i]; }
406
407 Vector xPos_;
408 Vector yPos_;
409 BlockVector m_;
410};
411
412} // end namespace Dumux
413
414#endif
The common code for all 3rd order polynomial splines.
Definition: adapt.hh:29
The common code for all 3rd order polynomial splines with more than two sampling points.
Definition: fixedlengthspline_.hh:44
FixedLengthSpline_()
Definition: fixedlengthspline_.hh:53
BlockVector m_
Definition: fixedlengthspline_.hh:409
void setContainerOfTuples(const XYContainer &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition: fixedlengthspline_.hh:197
void makeFullSpline_(Scalar m0, Scalar m1)
Create a full spline from the already set sampling points.
Definition: fixedlengthspline_.hh:359
void setXYArrays(int numberSamples, const ScalarArrayX &x, const ScalarArrayY &y)
Set the sampling points natural spline using C-style arrays.
Definition: fixedlengthspline_.hh:232
Vector xPos_
Definition: fixedlengthspline_.hh:407
void setContainerOfPoints(const XYContainer &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition: fixedlengthspline_.hh:164
void makeNaturalSpline_()
Create a natural spline from the already set sampling points.
Definition: fixedlengthspline_.hh:376
void setXYArrays(int numberSamples, const ScalarArrayX &x, const ScalarArrayY &y, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using C-style arrays.
Definition: fixedlengthspline_.hh:84
void setArrayOfPoints(int numberSamples, const PointArray &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using a C-style array.
Definition: fixedlengthspline_.hh:135
Vector yPos_
Definition: fixedlengthspline_.hh:408
Scalar x_(int i) const
Returns the x coordinate of the i-th sampling point.
Definition: fixedlengthspline_.hh:391
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y)
Set the sampling points of a natural spline using STL-compatible containers.
Definition: fixedlengthspline_.hh:258
void setContainerOfPoints(const XYContainer &points)
Set the sampling points of a natural spline using a STL-compatible container of array-like objects.
Definition: fixedlengthspline_.hh:308
Scalar y_(int i) const
Returns the y coordinate of the i-th sampling point.
Definition: fixedlengthspline_.hh:397
Scalar moment_(int i) const
Returns the moment (i.e. second derivative) of the spline at the i-th sampling point.
Definition: fixedlengthspline_.hh:404
int numSamples() const
Returns the number of sampling points.
Definition: fixedlengthspline_.hh:61
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using STL-compatible containers.
Definition: fixedlengthspline_.hh:110
void setArrayOfPoints(int numberSamples, const PointArray &points)
Set the sampling points of a natural spline using a C-style array.
Definition: fixedlengthspline_.hh:284
void setContainerOfTuples(const XYContainer &points)
Set the sampling points of a natural spline using a STL-compatible container of tuple-like objects.
Definition: fixedlengthspline_.hh:338
The common code for all 3rd order polynomial splines.
Definition: splinecommon_.hh:44
void makeFullSystem_(Matrix &M, Vector &d, Scalar m0, Scalar m1)
Make the linear system of equations Mx = d which results in the moments of the full spline.
Definition: splinecommon_.hh:427
void makeNaturalSystem_(Matrix &M, Vector &d)
Make the linear system of equations Mx = d which results in the moments of the natural spline....
Definition: splinecommon_.hh:452