From 3fcf5296a636a642620661293e090466d84ea1d6 Mon Sep 17 00:00:00 2001 From: Andreas Lauser <and@poware.org> Date: Thu, 16 Jun 2011 15:53:05 +0000 Subject: [PATCH] spline: rework the set() methods the old set() methods are still around but deprecated. Everything should still compile... git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@6000 2fb0f335-1f38-0410-981e-8018bf24f1b0 --- dumux/common/fixedlengthspline_.hh | 364 +++++++++++++++++++----- dumux/common/spline.hh | 234 +++++++++++----- dumux/common/variablelengthspline_.hh | 384 +++++++++++++++++++------- test/common/spline/test_spline.cc | 38 +-- 4 files changed, 763 insertions(+), 257 deletions(-) diff --git a/dumux/common/fixedlengthspline_.hh b/dumux/common/fixedlengthspline_.hh index 8fb08a47ae..8275bfc224 100644 --- a/dumux/common/fixedlengthspline_.hh +++ b/dumux/common/fixedlengthspline_.hh @@ -60,129 +60,343 @@ public: int numSamples() const { return nSamples; } + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// + // Full splines // + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// + /*! - * \brief Set the sampling points and the boundary slopes of the - * spline. + * \brief Set the sampling points and the boundary slopes of a + * full spline using C-style arrays. + * + * This method uses separate array-like objects for the values of + * the X and Y coordinates. In this context 'array-like' means + * that an access to the members is provided via the [] + * operator. (e.g. C arrays, std::vector, std::array, etc.) Each + * array must be of size 'numberSamples' at least. Also, the number of + * sampling points must be larger than 1. */ + template <class ScalarArrayX, class ScalarArrayY> + void setXYArrays(int numberSamples, + const ScalarArrayX &x, + const ScalarArrayY &y, + Scalar m0, Scalar m1) + { + assert(numberSamples == numSamples()); + + for (int i = 0; i < numberSamples; ++i) { + xPos_[i] = x[i]; + yPos_[i] = y[i]; + } + makeFullSpline_(m0, m1); + } + template <class ScalarContainer> + DUNE_DEPRECATED // use setXYArrays void set(const ScalarContainer &x, const ScalarContainer &y, - Scalar m0, Scalar m1) + Scalar m0, + Scalar m1) + { setXYArrays(numSamples(), x, y, m0, m1); } + + template <class ScalarArray> + DUNE_DEPRECATED // use setXYArrays + void set(int numberSamples, + const ScalarArray &x, + const ScalarArray &y, + Scalar m0, + Scalar m1) + { setXYArrays(numberSamples, x, y, m0, m1); } + + /*! + * \brief Set the sampling points and the boundary slopes of a + * full spline using STL-compatible containers. + * + * This method uses separate STL-compatible containers for the + * values of the sampling points' X and Y + * coordinates. "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. + */ + template <class ScalarContainerX, class ScalarContainerY> + void setXYContainers(const ScalarContainerX &x, + const ScalarContainerY &y, + Scalar m0, Scalar m1) + { + assert(x.size() == y.size()); + assert(x.size() == numSamples()); + + std::copy(x.begin(), x.end(), xPos_.begin()); + std::copy(y.begin(), y.end(), yPos_.begin()); + makeFullSpline_(m0, m1); + } + + /*! + * \brief Set the sampling points and the boundary slopes of a + * full spline using a C-style array. + * + * This method uses a single array of sampling points, which are + * seen as an array-like object which provides access to the X and + * Y coordinates. In this context 'array-like' means that an + * access to the members is provided via the [] operator. (e.g. C + * arrays, std::vector, std::array, etc.) The array containing + * the sampling points must be of size 'numberSamples' at least. Also, + * the number of sampling points must be larger than 1. + */ + template <class PointArray> + void setArrayOfPoints(int numberSamples, + const PointArray &points, + Scalar m0, + Scalar m1) { - this->assignSamplingPoints_(xPos_, yPos_, x, y, numSamples()); - myMakeFullSpline_(m0, m1); + // a spline with no or just one sampling points? what an + // incredible bad idea! + assert(numberSamples == numSamples()); + + for (int i = 0; i < numberSamples; ++i) { + xPos_[i] = points[i][0]; + yPos_[i] = points[i][1]; + } + makeFullSpline_(m0, m1); } + template <class XYContainer> + DUNE_DEPRECATED // use setArrayOfPoints instead + void set(const XYContainer &points, + Scalar m0, + Scalar m1) + { setArrayOfPoints(numSamples(), points, m0, m1); } + + /*! - * \brief Set the sampling points using an array of (x,y) pairs. + * \brief Set the sampling points and the boundary slopes of a + * full spline using a STL-compatible container of + * array-like objects. * - * Variant for a full spline. + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be array-like objects storing the + * X and Y coordinates. "STL-compatible" means that the container + * provides access to iterators using the begin(), end() methods + * and also provides a size() method. Also, the number of entries + * in the X and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaArrayContainer(const XYContainer &points, Scalar m0, Scalar m1) + void setContainerOfPoints(const XYContainer &points, + Scalar m0, + Scalar m1) { assert(points.size() == numSamples()); - this->assignFromArrayList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); - myMakeFullSpline_(m0, m1); + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++i, ++it) { + xPos_[i] = (*it)[0]; + yPos_[i] = (*it)[1]; + } + + // make a full spline + makeFullSpline_(m0, m1); } /*! - * \brief Set the sampling points using an array of (x,y) pairs. + * \brief Set the sampling points and the boundary slopes of a + * full spline using a STL-compatible container of + * tuple-like objects. * - * Variant for a full spline. + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be tuple-like objects storing the + * X and Y coordinates. "tuple-like" means that the objects + * provide access to the x values via std::get<0>(obj) and to the + * y value via std::get<1>(obj) (e.g. std::tuple or + * std::pair). "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaArrayArray(int nuSamples, const XYContainer &points, Scalar m0, Scalar m1) + void setContainerOfTuples(const XYContainer &points, + Scalar m0, + Scalar m1) { - assert(nuSamples == numSamples()); - - this->assignFromArrayList_(xPos_, - yPos_, - &points[0], - &points[numSamples() - 1], - numSamples()); - myMakeFullSpline_(m0, m1); - } + assert(points.size() == numSamples()); + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++i, ++it) { + xPos_[i] = std::get<0>(*it); + yPos_[i] = std::get<1>(*it); + } + + // make a full spline + makeFullSpline_(m0, m1); + } + + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// + // Natural splines // + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// /*! - * \brief Set the sampling points using an array of (x,y) tuples. + * \brief Set the sampling points natural spline using C-style arrays. * - * Variant for a full spline. + * This method uses separate array-like objects for the values of + * the X and Y coordinates. In this context 'array-like' means + * that an access to the members is provided via the [] + * operator. (e.g. C arrays, std::vector, std::array, etc.) Each + * array must be of size 'numberSamples' at least. Also, the number of + * sampling points must be larger than 1. */ - template <class TupleContainer> - void setViaTupleContainer(const TupleContainer &points, Scalar m0, Scalar m1) + template <class ScalarArrayX, class ScalarArrayY> + void setXYArrays(int numberSamples, + const ScalarArrayX &x, + const ScalarArrayY &y) { - this->assignFromTupleList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); - myMakeFullSpline_(m0, m1); + assert(numberSamples == numSamples()); + + for (int i = 0; i < numberSamples; ++i) { + xPos_[i] = x[i]; + yPos_[i] = y[i]; + } + + makeNaturalSpline_(); } + template <class ScalarArray> + DUNE_DEPRECATED // use setXYArrays + void set(int numberSamples, + const ScalarArray &x, + const ScalarArray &y) + { setXYArrays(numberSamples, x, y); } + + template <class ScalarContainer> + DUNE_DEPRECATED // use setXYArrays + void set(const ScalarContainer &x, + const ScalarContainer &y) + { setXYArrays(numSamples(), x, y); } + /*! - * \brief Set the sampling points for a natural spline. + * \brief Set the sampling points of a natural spline using + * STL-compatible containers. + * + * This method uses separate STL-compatible containers for the + * values of the sampling points' X and Y + * coordinates. "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ - template <class ScalarContainer> - void set(const ScalarContainer &x, const ScalarContainer &y) + template <class ScalarContainerX, class ScalarContainerY> + void setXYContainers(const ScalarContainerX &x, + const ScalarContainerY &y) { - this->assignSamplingPoints_(xPos_, yPos_, x, y, numSamples()); - myMakeNaturalSpline_(); - } + assert(x.size() == y.size()); + assert(x.size() > 1); + setNumSamples_(x.size()); + std::copy(x.begin(), x.end(), xPos_.begin()); + std::copy(y.begin(), y.end(), yPos_.begin()); + makeNaturalSpline_(); + } + + /*! - * \brief Set the sampling points using an array of (x,y) tuples. + * \brief Set the sampling points of a natural spline using a + * C-style array. * - * Variant for a natural spline. + * This method uses a single array of sampling points, which are + * seen as an array-like object which provides access to the X and + * Y coordinates. In this context 'array-like' means that an + * access to the members is provided via the [] operator. (e.g. C + * arrays, std::vector, std::array, etc.) The array containing + * the sampling points must be of size 'numberSamples' at least. Also, + * the number of sampling points must be larger than 1. */ - template <class XYContainer> - void setViaArrayContainer(const XYContainer &points) + template <class PointArray> + void setArrayOfPoints(int numberSamples, + const PointArray &points) { - this->assignFromArrayList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); - myMakeNaturalSpline_(); + assert(numberSamples == numSamples()); + + for (int i = 0; i < numberSamples; ++i) { + xPos_[i] = points[i][0]; + yPos_[i] = points[i][1]; + } + makeNaturalSpline_(); } + template <class XYContainer> + DUNE_DEPRECATED // use setArrayOfPoints instead + void set(int numberSamples, + const XYContainer &points) + { setArrayOfPoints(numberSamples, points); } + + template <class XYArray> + DUNE_DEPRECATED // use setArrayOfPoints() instead + void set(const XYArray &points) + { setArrayOfPoints(numSamples(), points); } + + /*! - * \brief Set the sampling points using an array of (x,y) tuples. + * \brief Set the sampling points of a natural spline using a + * STL-compatible container of array-like objects. * - * Variant for a natural spline. + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be array-like objects storing the + * X and Y coordinates. "STL-compatible" means that the container + * provides access to iterators using the begin(), end() methods + * and also provides a size() method. Also, the number of entries + * in the X and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaArrayArray(int nuSamples, const XYContainer &points) + void setContainerOfPoints(const XYContainer &points) { - assert(nuSamples == numSamples()); - - this->assignFromArrayList_(xPos_, - yPos_, - &points[0], - &points[numSamples() - 1], - numSamples()); - myMakeNaturalSpline_(); + assert(points.size() == numSamples()); + + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++ i, ++it) { + xPos_[i] = (*it)[0]; + yPos_[i] = (*it)[1]; + } + + // make a natural spline + makeNaturalSpline_(); } /*! - * \brief Set the sampling points using an array of (x,y) tuples. + * \brief Set the sampling points of a natural spline using a + * STL-compatible container of tuple-like objects. * - * Variant for a natural spline. + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be tuple-like objects storing the + * X and Y coordinates. "tuple-like" means that the objects + * provide access to the x values via std::get<0>(obj) and to the + * y value via std::get<1>(obj) (e.g. std::tuple or + * std::pair). "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ - template <class TupleContainer> - void setViaTupleContainer(const TupleContainer &points) + template <class XYContainer> + void setContainerOfTuples(const XYContainer &points) { - this->assignFromTupleList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); - myMakeNaturalSpline_(); + assert(points.size() == numSamples()); + + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++i, ++it) { + xPos_[i] = std::get<0>(*it); + yPos_[i] = std::get<1>(*it); + } + + // make a natural spline + makeNaturalSpline_(); } protected: @@ -191,7 +405,7 @@ protected: * * Also creates temporary matrix and right hand side vector. */ - void myMakeFullSpline_(Scalar m0, Scalar m1) + void makeFullSpline_(Scalar m0, Scalar m1) { BTDMatrix M(numSamples()); BlockVector d(numSamples()); @@ -208,7 +422,7 @@ protected: * * Also creates temporary matrix and right hand side vector. */ - void myMakeNaturalSpline_() + void makeNaturalSpline_() { BTDMatrix M(numSamples()); BlockVector d(numSamples()); diff --git a/dumux/common/spline.hh b/dumux/common/spline.hh index 63b6a75dde..fd0baa76db 100644 --- a/dumux/common/spline.hh +++ b/dumux/common/spline.hh @@ -49,8 +49,6 @@ namespace Dumux s''(x_n) & = 0 \f} */ - - template<class Scalar, int numSamples = 2> class Spline : public FixedLengthSpline_<Scalar, numSamples> { @@ -64,23 +62,25 @@ public: { }; /*! - * \brief Convenience constructor for a natural spline + * \brief Convenience constructor for a full spline * * \param x An array containing the \f$x\f$ values of the spline's sampling points * \param y An array containing the \f$y\f$ values of the spline's sampling points */ - template <class ScalarContainer> - Spline(const ScalarContainer &x, - const ScalarContainer &y) - { this->set(x, y); } + template <class ScalarArray> + Spline(const ScalarArray &x, + const ScalarArray &y) + { this->setXYArrays(numSamples, x, y); } /*! - * \brief Convenience constructor for a natural spline + * \brief Convenience constructor for a full spline * - * \param tuples An array of \f$(x,y)\f$ tuples of the spline's sampling points + * \param x An array containing the \f$x\f$ values of the spline's sampling points + * \param y An array containing the \f$y\f$ values of the spline's sampling points */ - Spline(const std::initializer_list<const std::pair<Scalar, Scalar> > &points) - { this->setViaTupleContainer(points); } + template <class PointArray> + Spline(const PointArray &points) + { this->setArrayOfPoints(numSamples, points); } /*! * \brief Convenience constructor for a full spline @@ -90,12 +90,12 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class ScalarContainer> - Spline(const ScalarContainer &x, - const ScalarContainer &y, + template <class ScalarArray> + Spline(const ScalarArray &x, + const ScalarArray &y, Scalar m0, Scalar m1) - { this->set(x, y, m0, m1); } + { this->setXYArrays(numSamples, x, y, m0, m1); } /*! * \brief Convenience constructor for a full spline @@ -104,10 +104,11 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - Spline(const std::initializer_list<const std::pair<Scalar, Scalar> > &points, + template <class PointArray> + Spline(const PointArray &points, Scalar m0, Scalar m1) - { this->setViaTupleContainer(points, m0, m1); } + { this->setArrayOfPoints(numSamples, points, m0, m1); } }; /*! @@ -131,7 +132,7 @@ public: \f} */ template<class Scalar> -class Spline<Scalar, -1> : public VariableLengthSpline_<Scalar> +class Spline<Scalar, /*numSamples=*/-1> : public VariableLengthSpline_<Scalar> { public: /*! @@ -149,11 +150,11 @@ public: * \param x An array containing the \f$x\f$ values of the spline's sampling points * \param y An array containing the \f$y\f$ values of the spline's sampling points */ - template <class ScalarContainer> + template <class ScalarArrayX, class ScalarArrayY> Spline(int nSamples, - const ScalarContainer &x, - const ScalarContainer &y) - { this->set(nSamples, x, y); } + const ScalarArrayX &x, + const ScalarArrayY &y) + { this->setXYArrays(nSamples, x, y); } /*! * \brief Convenience constructor for a natural spline @@ -161,10 +162,10 @@ public: * \param nSamples The number of sampling points (must be > 2) * \param tuples An array of \f$(x,y)\f$ tuples of the spline's sampling points */ - template <class XYContainer> + template <class PointArray> Spline(int nSamples, - const XYContainer &tuples) - { this->set(nSamples, tuples); } + const PointArray &points) + { this->setArrayOfPoints(nSamples, points); } /*! * \brief Convenience constructor for a natural spline @@ -175,15 +176,18 @@ public: template <class ScalarContainer> Spline(const ScalarContainer &x, const ScalarContainer &y) - { this->set(x, y); } + { this->setXYContainers(x, y); } /*! * \brief Convenience constructor for a natural spline * - * \param tuples An array of \f$(x,y)\f$ tuples of the spline's sampling points (must have a size() method) + * \param points An array of \f$(x,y)\f$ tuples of the spline's sampling points (must have a size() method) + * \param m0 The slope of the spline at \f$x_0\f$ + * \param m1 The slope of the spline at \f$x_n\f$ */ - Spline(const std::initializer_list<const std::pair<Scalar, Scalar> > &points) - { this->setViaArrayContainer(points); } + template <class PointContainer> + Spline(const PointContainer &points) + { this->setContainerOfPoints(points); } /*! * \brief Convenience constructor for a full spline @@ -194,28 +198,28 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class ScalarContainer> + template <class ScalarArray> Spline(int nSamples, - const ScalarContainer &x, - const ScalarContainer &y, + const ScalarArray &x, + const ScalarArray &y, Scalar m0, Scalar m1) - { this->set(nSamples, x, y, m0, m1); } + { this->setXYArrays(nSamples, x, y, m0, m1); } /*! * \brief Convenience constructor for a full spline * * \param nSamples The number of sampling points (must be >= 2) - * \param points An array of \f$(x,y)\f$ tuples of the spline's sampling points + * \param points An array containing the \f$x\f$ and \f$x\f$ values of the spline's sampling points * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class XYContainer> + template <class PointArray> Spline(int nSamples, - const XYContainer &points, + const PointArray &points, Scalar m0, Scalar m1) - { this->setViaArrayArray(nSamples, points, m0, m1); } + { this->setArrayOfPoints(nSamples, points, m0, m1); } /*! * \brief Convenience constructor for a full spline @@ -225,12 +229,12 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class ScalarContainer> - Spline(const ScalarContainer &x, - const ScalarContainer &y, + template <class ScalarContainerX, class ScalarContainerY> + Spline(const ScalarContainerX &x, + const ScalarContainerY &y, Scalar m0, Scalar m1) - { this->set(x, y, m0, m1); } + { this->setXYContainers(x, y, m0, m1); } /*! * \brief Convenience constructor for a full spline @@ -239,18 +243,18 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class XYContainer> - Spline(const XYContainer &points, + template <class PointContainer> + Spline(const PointContainer &points, Scalar m0, Scalar m1) - { this->setViaTupleContainer(points, m0, m1); } + { this->setContainerOfPoints(points, m0, m1); } }; /*! * \brief Do not allow splines with zero sampling points */ template<class Scalar> -class Spline<Scalar, 0> +class Spline<Scalar, /*numSamples=*/0> // Splines with zero sampling points do not make sense! { private: Spline() { }; }; @@ -258,7 +262,7 @@ class Spline<Scalar, 0> * \brief Do not allow splines with one sampling point */ template<class Scalar> -class Spline<Scalar, 1> +class Spline<Scalar, /*numSamples=*/1> // Splines with one sampling point do not make sense! { private: Spline() { }; }; @@ -286,13 +290,11 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class ScalarContainer> - Spline(const ScalarContainer &x, - const ScalarContainer &y, + template <class ScalarArrayX, class ScalarArrayY> + Spline(const ScalarArrayX &x, + const ScalarArrayY &y, Scalar m0, Scalar m1) - { - set(x,y,m0,m1); - } + { setXYArrays(2, x, y, m0, m1); } /*! * \brief Convenience constructor for a full spline @@ -301,11 +303,11 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_n\f$ */ - template <class XYContainer> - Spline(const XYContainer &points, + template <class PointArray> + Spline(const PointArray &points, Scalar m0, Scalar m1) - { this->set(points, m0, m1); } + { this->setArrayOfPoints(2, points, m0, m1); } /*! * \brief Convenience constructor for a full spline @@ -366,16 +368,50 @@ public: * \param m1 The slope of the spline at \f$x_1\f$ */ template <class ScalarContainer> - void set(const ScalarContainer &x, const ScalarContainer &y, - Scalar m0, Scalar m1) + void setXYArrays(int nSamples, + const ScalarContainer &x, + const ScalarContainer &y, + Scalar m0, Scalar m1) { + assert(nSamples == 2); + set(x[0], x[1], y[0], y[1], m0, m1); + } + + template <class ScalarArray> + DUNE_DEPRECATED // use setXYArrays + void set(const ScalarArray &x, + const ScalarArray &y, + Scalar m0, + Scalar m1) + { setXYArrays(2, x, y, m0, m1); } + + /*! + * \brief Set the sampling points and the boundary slopes of the + * spline. + * + * \param x An array containing the \f$x\f$ values of the sampling points + * \param y An array containing the \f$y\f$ values of the sampling points + * \param m0 The slope of the spline at \f$x_0\f$ + * \param m1 The slope of the spline at \f$x_1\f$ + */ + template <class ScalarContainerX, class ScalarContainerY> + void setXYContainers(const ScalarContainerX &x, + const ScalarContainerY &y, + Scalar m0, Scalar m1) + { + assert(x.size() == y.size()); + assert(x.size() == 2); + Matrix M(numSamples()); Vector d; - assignXY_(x[0], x[1], y[0], y[1]); - this->makeFullSystem_(M, d, m0, m1); - // solve for the moments - M.solve(m_, d); + typename ScalarContainerX::const_iterator xIt0 = x.begin(); + typename ScalarContainerX::const_iterator xIt1 = xIt0; + ++xIt1; + typename ScalarContainerY::const_iterator yIt0 = y.begin(); + typename ScalarContainerY::const_iterator yIt1 = yIt0; + ++yIt1; + set(*xIt0, *xIt1, *yIt0, *yIt1); } /*! @@ -386,16 +422,80 @@ public: * \param m0 The slope of the spline at \f$x_0\f$ * \param m1 The slope of the spline at \f$x_1\f$ */ - template <class XYContainer> - void set(const XYContainer &points, Scalar m0, Scalar m1) + template <class PointArray> + void setArrayOfPoints(int nSamples, + const PointArray &points, + Scalar m0, + Scalar m1) { + assert(nSamples == 2); + + set(points[0][0], + points[1][0], + points[0][1], + points[1][1], + m0, m1); + } + + template <class PointArray> + DUNE_DEPRECATED // use setArrayOfPoints + void set(const PointArray &points, + Scalar m0, + Scalar m1) + { setArrayOfPoints(2, points, m0, m1); } + + /*! + * \brief Set the sampling points and the boundary slopes from an + * STL-like container of points. + * + * \param points An array of \f$(x,y)\f$ tuples of the spline's sampling points + * \param m0 The slope of the spline at \f$x_0\f$ + * \param m1 The slope of the spline at \f$x_1\f$ + */ + template <class PointContainer> + void setContainerOfPoints(const PointContainer &points, + Scalar m0, + Scalar m1) + { + assert(points.size() == 2); + Matrix M; Vector d; - assignXY_(points[0][0], points[1][0], points[0][1], points[1][1]); - this->makeFullSystem_(M, d, m0, m1); + typename PointContainer::const_iterator it0 = points.begin(); + typename PointContainer::const_iterator it1 = it0; + ++it1; + + set((*it0)[0], + (*it0)[1], + (*it1)[0], + (*it1)[1], + m0, m1); + } - // solve for the moments - M.solve(m_, d); + /*! + * \brief Set the sampling points and the boundary slopes from an + * STL-like container of tuples. + * + * \param points An array of \f$(x,y)\f$ tuples of the spline's sampling points + * \param m0 The slope of the spline at \f$x_0\f$ + * \param m1 The slope of the spline at \f$x_1\f$ + */ + template <class TupleContainer> + void setContainerOfTuples(const TupleContainer &tuples, + Scalar m0, + Scalar m1) + { + assert(tuples.size() == 2); + + typename TupleContainer::const_iterator it0 = tuples.begin(); + typename TupleContainer::const_iterator it1 = it0; + ++it1; + + set(std::get<0>(*it0), + std::get<1>(*it0), + std::get<0>(*it1), + std::get<1>(*it1), + m0, m1); } protected: diff --git a/dumux/common/variablelengthspline_.hh b/dumux/common/variablelengthspline_.hh index 1537618aa8..592b24b1ef 100644 --- a/dumux/common/variablelengthspline_.hh +++ b/dumux/common/variablelengthspline_.hh @@ -55,176 +55,368 @@ public: int numSamples() const { return xPos_.size(); } + + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// + // Full splines // + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// + /*! - * \brief Set the sampling points and the boundary slopes of the - * spline. + * \brief Set the sampling points and the boundary slopes of a + * full spline using C-style arrays. + * + * This method uses separate array-like objects for the values of + * the X and Y coordinates. In this context 'array-like' means + * that an access to the members is provided via the [] + * operator. (e.g. C arrays, std::vector, std::array, etc.) Each + * array must be of size 'nSamples' at least. Also, the number of + * sampling points must be larger than 1. */ - template <class ScalarContainer> - void set(int nSamples, - const ScalarContainer &x, - const ScalarContainer &y, - Scalar m0, Scalar m1) + template <class ScalarArrayX, class ScalarArrayY> + void setXYArrays(int nSamples, + const ScalarArrayX &x, + const ScalarArrayY &y, + Scalar m0, Scalar m1) { + assert(nSamples > 1); + setNumSamples_(nSamples); - this->assignSamplingPoints_(xPos_, yPos_, x, y, numSamples()); - myMakeFullSpline_(m0, m1); + for (int i = 0; i < nSamples; ++i) { + xPos_[i] = x[i]; + yPos_[i] = y[i]; + } + makeFullSpline_(m0, m1); } + template <class ScalarArray> + DUNE_DEPRECATED // use setXYArrays + void set(int nSamples, + const ScalarArray &x, + const ScalarArray &y, + Scalar m0, + Scalar m1) + { setXYArrays(nSamples, x, y, m0, m1); } + /*! - * \brief Set the sampling points and the boundary slopes of the - * spline. + * \brief Set the sampling points and the boundary slopes of a + * full spline using STL-compatible containers. * - * This assumes that the ScalarContainer type has a size method. + * This method uses separate STL-compatible containers for the + * values of the sampling points' X and Y + * coordinates. "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ - template <class ScalarContainer> - void set(const ScalarContainer &x, - const ScalarContainer &y, - Scalar m0, Scalar m1) + template <class ScalarContainerX, class ScalarContainerY> + void setXYContainers(const ScalarContainerX &x, + const ScalarContainerY &y, + Scalar m0, Scalar m1) { assert(x.size() == y.size()); - set(x.size(), x, y, m0, m1); + assert(x.size() > 1); + + setNumSamples_(x.size()); + std::copy(x.begin(), x.end(), xPos_.begin()); + std::copy(y.begin(), y.end(), yPos_.begin()); + makeFullSpline_(m0, m1); } + + template <class ScalarContainer> + DUNE_DEPRECATED // use setXYContainers + void set(const ScalarContainer &x, + const ScalarContainer &y, + Scalar m0, + Scalar m1) + { setXYContainers(x, y, m0, m1); } + /*! - * \brief Set the sampling points for a full spline using a - * container of (x,y) tuples. + * \brief Set the sampling points and the boundary slopes of a + * full spline using a C-style array. + * + * This method uses a single array of sampling points, which are + * seen as an array-like object which provides access to the X and + * Y coordinates. In this context 'array-like' means that an + * access to the members is provided via the [] operator. (e.g. C + * arrays, std::vector, std::array, etc.) The array containing + * the sampling points must be of size 'nSamples' at least. Also, + * the number of sampling points must be larger than 1. */ - template <class XYContainer> - void setViaArrayArray(int nSamples, - const XYContainer &points, + template <class PointArray> + void setArrayOfPoints(int nSamples, + const PointArray &points, Scalar m0, Scalar m1) { + // a spline with no or just one sampling points? what an + // incredible bad idea! + assert(nSamples > 1); + setNumSamples_(nSamples); - this->assignFromArrayList_(xPos_, - yPos_, - &points[0], - &points[numSamples() - 1], - numSamples()); - myMakeFullSpline_(m0, m1); + for (int i = 0; i < nSamples; ++i) { + xPos_[i] = points[i][0]; + yPos_[i] = points[i][1]; + } + makeFullSpline_(m0, m1); } + template <class XYContainer> + DUNE_DEPRECATED // use setArrayOfPoints instead + void set(int nSamples, + const XYContainer &points, + Scalar m0, + Scalar m1) + { setArrayOfPoints(nSamples, points, m0, m1); } + + /*! - * \brief Set the sampling points for a full spline using a - * container of (x,y) tuples. + * \brief Set the sampling points and the boundary slopes of a + * full spline using a STL-compatible container of + * array-like objects. * - * This assimes the XYContainer type has a size() method + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be array-like objects storing the + * X and Y coordinates. "STL-compatible" means that the container + * provides access to iterators using the begin(), end() methods + * and also provides a size() method. Also, the number of entries + * in the X and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaArrayContainer(const XYContainer &points, + void setContainerOfPoints(const XYContainer &points, Scalar m0, Scalar m1) - { + // a spline with no or just one sampling points? what an + // incredible bad idea! + assert(points.size() > 1); + setNumSamples_(points.size()); - this->assignFromArrayList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++i, ++it) { + xPos_[i] = (*it)[0]; + yPos_[i] = (*it)[1]; + } // make a full spline - myMakeFullSpline_(m0, m1); + makeFullSpline_(m0, m1); } + template <class XYContainer> + DUNE_DEPRECATED // use setContainerOfPoints() instead + void set(const XYContainer &points, + Scalar m0, + Scalar m1) + { setContainerOfPoints(points, m0, m1); } + /*! - * \brief Set the sampling points for a fill spline. + * \brief Set the sampling points and the boundary slopes of a + * full spline using a STL-compatible container of + * tuple-like objects. + * + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be tuple-like objects storing the + * X and Y coordinates. "tuple-like" means that the objects + * provide access to the x values via std::get<0>(obj) and to the + * y value via std::get<1>(obj) (e.g. std::tuple or + * std::pair). "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaTupleContainer(const XYContainer &points, + void setContainerOfTuples(const XYContainer &points, Scalar m0, Scalar m1) { // resize internal arrays setNumSamples_(points.size()); - - // assign the values to the internal arrays - this->assignFromTupleList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++i, ++it) { + xPos_[i] = std::get<0>(*it); + yPos_[i] = std::get<1>(*it); + } // make a full spline - myMakeFullSpline_(m0, m1); + makeFullSpline_(m0, m1); } - + + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// + // Natural splines // + /////////////////////////////////////// + /////////////////////////////////////// + /////////////////////////////////////// /*! - * \brief Set the sampling points for a natural spline. + * \brief Set the sampling points natural spline using C-style arrays. + * + * This method uses separate array-like objects for the values of + * the X and Y coordinates. In this context 'array-like' means + * that an access to the members is provided via the [] + * operator. (e.g. C arrays, std::vector, std::array, etc.) Each + * array must be of size 'nSamples' at least. Also, the number of + * sampling points must be larger than 1. */ - template <class ScalarContainer> - void set(int nSamples, - const ScalarContainer &x, - const ScalarContainer &y) + template <class ScalarArrayX, class ScalarArrayY> + void setXYArrays(int nSamples, + const ScalarArrayX &x, + const ScalarArrayY &y) { - // resize internal arrays + assert(nSamples > 1); + setNumSamples_(nSamples); - this->assignSamplingPoints_(xPos_, yPos_, x, y, numSamples()); - // make a natural spline - myMakeNaturalSpline_(); + for (int i = 0; i < nSamples; ++i) { + xPos_[i] = x[i]; + yPos_[i] = y[i]; + } + + makeNaturalSpline_(); } + template <class ScalarArray> + DUNE_DEPRECATED // use setXYArrays + void set(int nSamples, + const ScalarArray &x, + const ScalarArray &y) + { setXYArrays(nSamples, x, y); } + /*! - * \brief Set the sampling points for a natural spline. + * \brief Set the sampling points of a natural spline using + * STL-compatible containers. * - * This assimes the ScalarContainer type has a size() method + * This method uses separate STL-compatible containers for the + * values of the sampling points' X and Y + * coordinates. "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ - template <class ScalarContainer> - void set(const ScalarContainer &x, - const ScalarContainer &y) + template <class ScalarContainerX, class ScalarContainerY> + void setXYContainers(const ScalarContainerX &x, + const ScalarContainerY &y) { assert(x.size() == y.size()); - set(x.size(), x, y); + assert(x.size() > 1); + + setNumSamples_(x.size()); + std::copy(x.begin(), x.end(), xPos_.begin()); + std::copy(y.begin(), y.end(), yPos_.begin()); + makeNaturalSpline_(); } + + template <class ScalarContainer> + DUNE_DEPRECATED // use setXYContainers + void set(const ScalarContainer &x, + const ScalarContainer &y) + { setXYContainers(x, y); } + /*! - * \brief Set the sampling points for a natural spline using a container of (x,y) tuples. + * \brief Set the sampling points of a natural spline using a + * C-style array. + * + * This method uses a single array of sampling points, which are + * seen as an array-like object which provides access to the X and + * Y coordinates. In this context 'array-like' means that an + * access to the members is provided via the [] operator. (e.g. C + * arrays, std::vector, std::array, etc.) The array containing + * the sampling points must be of size 'nSamples' at least. Also, + * the number of sampling points must be larger than 1. */ - template <class XYContainer> - void setViaArrayArray(int nSamples, - const XYContainer &points) + template <class PointArray> + void setArrayOfPoints(int nSamples, + const PointArray &points) { + // a spline with no or just one sampling points? what an + // incredible bad idea! + assert(nSamples > 1); + setNumSamples_(nSamples); - this->assignFromArrayList_(xPos_, - yPos_, - &points[0], - &points[numSamples() - 1], - numSamples()); - myMakeNaturalSpline_(); + for (int i = 0; i < nSamples; ++i) { + xPos_[i] = points[i][0]; + yPos_[i] = points[i][1]; + } + makeNaturalSpline_(); } + template <class XYContainer> + DUNE_DEPRECATED // use setArrayOfPoints instead + void set(int nSamples, + const XYContainer &points) + { setArrayOfPoints(nSamples, points); } + + /*! - * \brief Set the sampling points for a natural spline using a container of (x,y) tuples. + * \brief Set the sampling points of a natural spline using a + * STL-compatible container of array-like objects. * - * This assimes the XYContainer type has a size() method + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be array-like objects storing the + * X and Y coordinates. "STL-compatible" means that the container + * provides access to iterators using the begin(), end() methods + * and also provides a size() method. Also, the number of entries + * in the X and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaArrayContainer(const XYContainer &points) + void setContainerOfPoints(const XYContainer &points) { + // a spline with no or just one sampling points? what an + // incredible bad idea! + assert(points.size() > 1); + setNumSamples_(points.size()); - this->assignFromArrayList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); - myMakeNaturalSpline_(); + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++ i, ++it) { + xPos_[i] = (*it)[0]; + yPos_[i] = (*it)[1]; + } + + // make a natural spline + makeNaturalSpline_(); } + template <class XYContainer> + DUNE_DEPRECATED // use setContainerOfPoints() instead + void set(const XYContainer &points) + { setContainerOfPoints(points); } + /*! - * \brief Set the sampling points for a natural spline. + * \brief Set the sampling points of a natural spline using a + * STL-compatible container of tuple-like objects. + * + * This method uses a single STL-compatible container of sampling + * points, which are assumed to be tuple-like objects storing the + * X and Y coordinates. "tuple-like" means that the objects + * provide access to the x values via std::get<0>(obj) and to the + * y value via std::get<1>(obj) (e.g. std::tuple or + * std::pair). "STL-compatible" means that the container provides + * access to iterators using the begin(), end() methods and also + * provides a size() method. Also, the number of entries in the X + * and the Y containers must be equal and larger than 1. */ template <class XYContainer> - void setViaTupleContainer(const XYContainer &points) + void setContainerOfTuples(const XYContainer &points) { + // resize internal arrays setNumSamples_(points.size()); - this->assignFromTupleList_(xPos_, - yPos_, - points.begin(), - points.end(), - numSamples()); - myMakeNaturalSpline_(); - } + typename XYContainer::const_iterator it = points.begin(); + typename XYContainer::const_iterator endIt = points.end(); + for (int i = 0; it != endIt; ++i, ++it) { + xPos_[i] = std::get<0>(*it); + yPos_[i] = std::get<1>(*it); + } + // make a natural spline + makeNaturalSpline_(); + } + protected: /*! * \brief Resizes the internal vectors to store the sample points. @@ -237,11 +429,11 @@ protected: } /*! - * \brief Create a full spline from the already set sampling points. + * \brief Create a natural spline from the already set sampling points. * * Also creates temporary matrix and right hand side vector. */ - void myMakeFullSpline_(Scalar m0, Scalar m1) + void makeFullSpline_(Scalar m0, Scalar m1) { BTDMatrix M(numSamples()); Vector d(numSamples()); @@ -258,7 +450,7 @@ protected: * * Also creates temporary matrix and right hand side vector. */ - void myMakeNaturalSpline_() + void makeNaturalSpline_() { BTDMatrix M(numSamples()); Vector d(numSamples()); diff --git a/test/common/spline/test_spline.cc b/test/common/spline/test_spline.cc index c478b82316..2cd9fd9d05 100644 --- a/test/common/spline/test_spline.cc +++ b/test/common/spline/test_spline.cc @@ -157,8 +157,8 @@ void testAll() // full spline { Dumux::Spline<double, 2> sp(x[0], x[1], y[0], y[1], m0, m1); sp.set(x[0],x[1],y[0],y[1],m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, 2> sp(x, y, m0, m1); sp.set(x,y,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, 2> sp(points, m0, m1); sp.set(points,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, 2> sp(x, y, m0, m1); sp.setXYArrays(2, x, y, m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, 2> sp(points, m0, m1); sp.setArrayOfPoints(2, points, m0, m1); testFull(sp, x, y, m0, m1); }; ///////// @@ -166,33 +166,33 @@ void testAll() ///////// // full spline - { Dumux::Spline<double, 5> sp(x, y, m0, m1); sp.set(x,y,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, 5> sp; sp.setViaArrayArray(5, points, m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, 5> sp; sp.setViaArrayContainer(pointVec,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, 5> sp(pointsInitList, m0, m1); sp.setViaTupleContainer(pointsInitList, m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, 5> sp(x, y, m0, m1); sp.setXYArrays(5, x, y, m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, 5> sp; sp.setArrayOfPoints(5, points, m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, 5> sp; sp.setContainerOfPoints(pointVec,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, 5> sp; sp.setContainerOfTuples(pointsInitList, m0, m1); testFull(sp, x, y, m0, m1); }; // natural spline - { Dumux::Spline<double, 5> sp(x, y); sp.set(x, y); testNatural(sp, x, y); }; - { Dumux::Spline<double, 5> sp; sp.setViaArrayArray(5, pointVec); testNatural(sp, x, y); }; - { Dumux::Spline<double, 5> sp(pointsInitList); sp.setViaTupleContainer(pointsInitList); testNatural(sp, x, y); }; + { Dumux::Spline<double, 5> sp(x, y); sp.setXYArrays(5, x, y); testNatural(sp, x, y); }; + { Dumux::Spline<double, 5> sp; sp.setContainerOfPoints(pointVec); testNatural(sp, x, y); }; + { Dumux::Spline<double, 5> sp; sp.setContainerOfTuples(pointsInitList); testNatural(sp, x, y); }; ///////// // test variable length splines ///////// // full spline - { Dumux::Spline<double, -1> sp(5, x, y, m0, m1); sp.set(5,x,y,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, -1> sp(xVec, yVec, m0, m1); sp.set(xVec,yVec,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, -1> sp; sp.setViaArrayArray(5,points,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, -1> sp; sp.setViaArrayContainer(pointVec,m0, m1); testFull(sp, x, y, m0, m1); }; - { Dumux::Spline<double, -1> sp(pointsInitList, m0, m1); sp.setViaTupleContainer(pointsInitList,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, -1> sp(5, x, y, m0, m1); sp.setXYArrays(5,x,y,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, -1> sp(xVec, yVec, m0, m1); sp.setXYContainers(xVec,yVec,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, -1> sp; sp.setArrayOfPoints(5,points,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, -1> sp; sp.setContainerOfPoints(pointVec,m0, m1); testFull(sp, x, y, m0, m1); }; + { Dumux::Spline<double, -1> sp; sp.setContainerOfTuples(pointsInitList,m0, m1); testFull(sp, x, y, m0, m1); }; // natural spline - { Dumux::Spline<double, -1> sp(5, x, y); sp.set(5,x,y); testNatural(sp, x, y); }; - { Dumux::Spline<double, -1> sp(xVec, yVec); sp.set(xVec,yVec); testNatural(sp, x, y); }; - { Dumux::Spline<double, -1> sp; sp.setViaArrayArray(5,points); testNatural(sp, x, y); }; - { Dumux::Spline<double, -1> sp; sp.setViaArrayContainer(pointVec); testNatural(sp, x, y); }; - { Dumux::Spline<double, -1> sp(pointsInitList, m0, m1); sp.setViaTupleContainer(pointsInitList); testNatural(sp, x, y); }; + { Dumux::Spline<double, -1> sp(5, x, y); sp.setXYArrays(5,x,y); testNatural(sp, x, y); }; + { Dumux::Spline<double, -1> sp(xVec, yVec); sp.setXYContainers(xVec,yVec); testNatural(sp, x, y); }; + { Dumux::Spline<double, -1> sp; sp.setArrayOfPoints(5,points); testNatural(sp, x, y); }; + { Dumux::Spline<double, -1> sp; sp.setContainerOfPoints(pointVec); testNatural(sp, x, y); }; + { Dumux::Spline<double, -1> sp; sp.setContainerOfTuples(pointsInitList); testNatural(sp, x, y); }; } void plot() -- GitLab