ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
cubic_spline.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2021 by the ryujin authors
4//
5
6#pragma once
7
8#ifdef DEAL_II_WITH_GSL
9#include <gsl/gsl_spline.h>
10
11namespace ryujin
12{
29 {
30 public:
37 CubicSpline(const std::vector<double> &x,
38 const std::vector<double> &y) noexcept
39 : x_(x)
40 , y_(y)
41 {
42 AssertNothrow(x_.size() == y_.size(), dealii::ExcInternalError());
43 AssertNothrow(x_.size() >= 2, dealii::ExcInternalError());
44 AssertNothrow(std::is_sorted(x_.begin(), x_.end()),
45 dealii::ExcInternalError());
46
47 spline = gsl_spline_alloc(gsl_interp_cspline, x_.size());
48 gsl_spline_init(spline, x_.data(), y_.data(), x_.size());
49
50 accel = gsl_interp_accel_alloc();
51 }
52
57 : CubicSpline(copy.x_, copy.y_)
58 {
59 }
60
64 CubicSpline &operator=(const CubicSpline &) = delete;
65
70 {
71 gsl_interp_accel_free(accel);
72 gsl_spline_free(spline);
73 }
74
81 inline double eval(double x) const
82 {
83 return gsl_spline_eval(spline, x, accel);
84 }
85
86 private:
87 const std::vector<double> x_;
88 const std::vector<double> y_;
89 gsl_spline *spline;
90 mutable gsl_interp_accel *accel;
91 };
92} // namespace ryujin
93
94#endif
CubicSpline & operator=(const CubicSpline &)=delete
CubicSpline(const CubicSpline &copy)
Definition: cubic_spline.h:56
double eval(double x) const
Definition: cubic_spline.h:81
CubicSpline(const std::vector< double > &x, const std::vector< double > &y) noexcept
Definition: cubic_spline.h:37