Line data Source code
1 0 : // Distributed under the MIT License. 2 : // See LICENSE.txt for details. 3 : 4 : #pragma once 5 : 6 : #include <gsl/gsl_spline.h> 7 : #include <memory> 8 : #include <vector> 9 : 10 : /// \cond 11 : namespace PUP { 12 : class er; 13 : } // namespace PUP 14 : /// \endcond 15 : 16 : namespace intrp { 17 : /*! 18 : * \ingroup NumericalAlgorithmsGroup 19 : * \brief A natural cubic spline interpolation class 20 : * 21 : * The class builds a cubic spline interpolant with natural boundary conditions 22 : * using the `x_values` and `y_values` passed into the constructor. For details 23 : * on the algorithm see the GSL documentation on `gsl_interp_cspline`. 24 : * 25 : * Here is an example how to use this class: 26 : * 27 : * \snippet Test_CubicSpline.cpp interpolate_example 28 : */ 29 1 : class CubicSpline { 30 : public: 31 0 : CubicSpline(std::vector<double> x_values, std::vector<double> y_values); 32 : 33 0 : CubicSpline() = default; 34 0 : CubicSpline(const CubicSpline& /*rhs*/); 35 0 : CubicSpline& operator=(const CubicSpline& /*rhs*/); 36 0 : CubicSpline(CubicSpline&& /*rhs*/) = default; 37 0 : CubicSpline& operator=(CubicSpline&& rhs) = default; 38 0 : ~CubicSpline() = default; 39 : 40 0 : bool operator==(const CubicSpline& rhs) const; 41 : 42 0 : double operator()(double x_to_interp_to) const; 43 : 44 0 : const std::vector<double>& x_values() const { return x_values_; } 45 : 46 0 : const std::vector<double>& y_values() const { return y_values_; } 47 : 48 : // NOLINTNEXTLINE(google-runtime-references) 49 0 : void pup(PUP::er& p); 50 : 51 : private: 52 0 : struct gsl_interp_accel_deleter { 53 0 : void operator()(gsl_interp_accel* acc) const; 54 : }; 55 0 : struct gsl_spline_deleter { 56 0 : void operator()(gsl_spline* spline) const; 57 : }; 58 : 59 0 : void initialize_interpolant(); 60 : 61 0 : std::vector<double> x_values_; 62 0 : std::vector<double> y_values_; 63 0 : std::unique_ptr<gsl_interp_accel, gsl_interp_accel_deleter> acc_; 64 0 : std::unique_ptr<gsl_spline, gsl_spline_deleter> spline_; 65 : }; 66 : } // namespace intrp