DFT-EFE
 
Loading...
Searching...
No Matches
Spline.h
Go to the documentation of this file.
1/*
2 * spline.h
3 *
4 * simple cubic spline interpolation library without external
5 * dependencies
6 *
7 * ---------------------------------------------------------------------
8 * Copyright (C) 2011, 2014, 2016, 2021 Tino Kluge (ttk448 at gmail.com)
9 *
10 * This program is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU General Public License
12 * as published by the Free Software Foundation; either version 2
13 * of the License, or (at your option) any later version.
14 *
15 * This program is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with this program. If not, see <http://www.gnu.org/licenses/>.
22 * ---------------------------------------------------------------------
23 *
24 */
25
26
27#ifndef dftefeSpline_h
28#define dftefeSpline_h
29
30#include <cstdio>
31#include <cassert>
32#include <cmath>
33#include <vector>
34#include <algorithm>
35#include <sstream>
36#include <string>
37
38
39// header file and we don't want to export symbols to the obj files
40namespace dftefe
41{
42 namespace utils
43 {
44 // spline interpolation
45 class Spline
46 {
47 public:
48 // spline types
50 {
51 linear = 10, // linear interpolation
52 cspline = 30, // cubic splines (classical C^2)
53 cspline_hermite = 31 // cubic hermite splines (local, only C^1)
54 };
55
56 // boundary condition type for the spline end-points
58 {
60 second_deriv = 2
61 };
62
63 private:
64 std::vector<double> d_x, d_y; // x,y coordinates of points
65 // interpolation parameters
66 // f(x) = a_i + b_i*(x-x_i) + c_i*(x-x_i)^2 + d_i*(x-x_i)^3
67 // where a_i = y_i, or else it won't go through grid points
68 std::vector<double> d_b, d_c, d_d; // spline coefficients
69 double d_c0; // for left extrapolation
75 double d_a, d_r;
76 unsigned int d_numSubDiv;
77 void
78 set_coeffs_from_b(); // calculate c_i, d_i from b_i
79 size_t
80 find_closest(double x) const; // closest idx so that d_x[idx]<=x
81
82 public:
83 // default constructor: set boundary condition to be zero curvature
84 // at both ends, i.e. natural splines
85 Spline();
86
87 Spline(const std::vector<double> &X,
88 const std::vector<double> &Y,
89 const bool isSubdivPowerLawGrid = false,
90 spline_type type = cspline,
91 bool make_monotonic = false,
92 bd_type left = second_deriv,
93 double left_value = 0.0,
94 bd_type right = second_deriv,
95 double right_value = 0.0);
96
97 // modify boundary conditions: if called it must be before set_points()
98 void
100 double left_value,
101 bd_type right,
102 double right_value);
103
104 // set all data points (cubic_spline=false means linear interpolation)
105 void
106 set_points(const std::vector<double> &x,
107 const std::vector<double> &y,
108 spline_type type = cspline);
109
110 // adjust coefficients so that the spline becomes piecewise monotonic
111 // where possible
112 // this is done by adjusting slopes at grid points by a non-negative
113 // factor and this will break C^2
114 // this can also break boundary conditions if adjustments need to
115 // be made at the boundary points
116 // returns false if no adjustments have been made, true otherwise
117 bool
119
120 // evaluates the spline at point x
121 double
122 operator()(double x) const;
123 std::vector<double>
124 coefficients(double x) const;
125 double
126 deriv(int order, double x) const;
127
128 // returns the input data points
129 std::vector<double>
130 get_x() const
131 {
132 return d_x;
133 }
134 std::vector<double>
135 get_y() const
136 {
137 return d_y;
138 }
139 double
140 get_x_min() const
141 {
142 assert(!d_x.empty());
143 return d_x.front();
144 }
145 double
146 get_x_max() const
147 {
148 assert(!d_x.empty());
149 return d_x.back();
150 }
151
152 // spline info string, i.e. spline type, boundary conditions etc.
153 std::string
154 info() const;
155 };
156
157 namespace splineInternal
158 {
159 // band matrix solver
161 {
162 private:
163 std::vector<std::vector<double>> d_upper; // upper band
164 std::vector<std::vector<double>> d_lower; // lower band
165 public:
166 band_matrix(){}; // constructor
167 band_matrix(int dim, int n_u, int n_l); // constructor
168 ~band_matrix(){}; // destructor
169 void
170 resize(int dim, int n_u, int n_l); // init with dim,n_u,n_l
171 int
172 dim() const; // matrix dimension
173 int
174 num_upper() const
175 {
176 return (int)d_upper.size() - 1;
177 }
178 int
179 num_lower() const
180 {
181 return (int)d_lower.size() - 1;
182 }
183 // access operator
184 double &
185 operator()(int i, int j); // write
186 double
187 operator()(int i, int j) const; // read
188 // we can store an additional diagonal (in d_lower)
189 double &
190 saved_diag(int i);
191 double
192 saved_diag(int i) const;
193 void
194 lu_decompose();
195 std::vector<double>
196 r_solve(const std::vector<double> &b) const;
197 std::vector<double>
198 l_solve(const std::vector<double> &b) const;
199 std::vector<double>
200 lu_solve(const std::vector<double> &b, bool is_lu_decomposed = false);
201 };
202 } // namespace splineInternal
203 } // namespace utils
204} // namespace dftefe
205
206#endif // dftefeSpline_h
Definition: Spline.h:46
double d_right_value
Definition: Spline.h:72
double d_left_value
Definition: Spline.h:72
double deriv(int order, double x) const
Definition: Spline.cpp:560
double get_x_max() const
Definition: Spline.h:146
spline_type d_type
Definition: Spline.h:70
void set_points(const std::vector< double > &x, const std::vector< double > &y, spline_type type=cspline)
Definition: Spline.cpp:192
double get_x_min() const
Definition: Spline.h:140
double d_r
Definition: Spline.h:75
unsigned int d_numSubDiv
Definition: Spline.h:76
double d_c0
Definition: Spline.h:69
void set_boundary(bd_type left, double left_value, bd_type right, double right_value)
Definition: Spline.cpp:152
std::vector< double > get_y() const
Definition: Spline.h:135
bd_type
Definition: Spline.h:58
@ second_deriv
Definition: Spline.h:60
@ first_deriv
Definition: Spline.h:59
bd_type d_left
Definition: Spline.h:71
spline_type
Definition: Spline.h:50
@ cspline
Definition: Spline.h:52
@ cspline_hermite
Definition: Spline.h:53
@ linear
Definition: Spline.h:51
std::vector< double > d_b
Definition: Spline.h:68
std::vector< double > d_y
Definition: Spline.h:64
double operator()(double x) const
Definition: Spline.cpp:490
bd_type d_right
Definition: Spline.h:71
std::vector< double > coefficients(double x) const
Definition: Spline.cpp:521
std::vector< double > d_d
Definition: Spline.h:68
bool d_made_monotonic
Definition: Spline.h:73
std::string info() const
Definition: Spline.cpp:623
size_t find_closest(double x) const
Definition: Spline.cpp:434
Spline()
Definition: Spline.cpp:111
std::vector< double > get_x() const
Definition: Spline.h:130
bool make_monotonic()
Definition: Spline.cpp:373
void set_coeffs_from_b()
Definition: Spline.cpp:166
std::vector< double > d_c
Definition: Spline.h:68
std::vector< double > d_x
Definition: Spline.h:64
bool d_isSubdivPowerLawGrid
Definition: Spline.h:74
double d_a
Definition: Spline.h:75
std::vector< std::vector< double > > d_lower
Definition: Spline.h:164
int num_lower() const
Definition: Spline.h:179
band_matrix()
Definition: Spline.h:166
std::vector< std::vector< double > > d_upper
Definition: Spline.h:163
~band_matrix()
Definition: Spline.h:168
double & saved_diag(int i)
Definition: Spline.cpp:711
std::vector< double > r_solve(const std::vector< double > &b) const
Definition: Spline.cpp:780
int num_upper() const
Definition: Spline.h:174
std::vector< double > l_solve(const std::vector< double > &b) const
Definition: Spline.cpp:762
int dim() const
Definition: Spline.cpp:664
void resize(int dim, int n_u, int n_l)
Definition: Spline.cpp:647
double & operator()(int i, int j)
Definition: Spline.cpp:680
std::vector< double > lu_solve(const std::vector< double > &b, bool is_lu_decomposed=false)
Definition: Spline.cpp:798
void lu_decompose()
Definition: Spline.cpp:719
dealii includes
Definition: AtomFieldDataSpherical.cpp:31