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
74 void
75 set_coeffs_from_b(); // calculate c_i, d_i from b_i
76 size_t
77 find_closest(double x) const; // closest idx so that d_x[idx]<=x
78
79 public:
80 // default constructor: set boundary condition to be zero curvature
81 // at both ends, i.e. natural splines
82 Spline();
83
84 Spline(const std::vector<double> &X,
85 const std::vector<double> &Y,
86 spline_type type = cspline,
87 bool make_monotonic = false,
88 bd_type left = second_deriv,
89 double left_value = 0.0,
90 bd_type right = second_deriv,
91 double right_value = 0.0);
92
93 // modify boundary conditions: if called it must be before set_points()
94 void
96 double left_value,
97 bd_type right,
98 double right_value);
99
100 // set all data points (cubic_spline=false means linear interpolation)
101 void
102 set_points(const std::vector<double> &x,
103 const std::vector<double> &y,
104 spline_type type = cspline);
105
106 // adjust coefficients so that the spline becomes piecewise monotonic
107 // where possible
108 // this is done by adjusting slopes at grid points by a non-negative
109 // factor and this will break C^2
110 // this can also break boundary conditions if adjustments need to
111 // be made at the boundary points
112 // returns false if no adjustments have been made, true otherwise
113 bool
115
116 // evaluates the spline at point x
117 double
118 operator()(double x) const;
119 std::vector<double>
120 coefficients(double x) const;
121 double
122 deriv(int order, double x) const;
123
124 // returns the input data points
125 std::vector<double>
126 get_x() const
127 {
128 return d_x;
129 }
130 std::vector<double>
131 get_y() const
132 {
133 return d_y;
134 }
135 double
136 get_x_min() const
137 {
138 assert(!d_x.empty());
139 return d_x.front();
140 }
141 double
142 get_x_max() const
143 {
144 assert(!d_x.empty());
145 return d_x.back();
146 }
147
148 // spline info string, i.e. spline type, boundary conditions etc.
149 std::string
150 info() const;
151 };
152
153 namespace splineInternal
154 {
155 // band matrix solver
157 {
158 private:
159 std::vector<std::vector<double>> d_upper; // upper band
160 std::vector<std::vector<double>> d_lower; // lower band
161 public:
162 band_matrix(){}; // constructor
163 band_matrix(int dim, int n_u, int n_l); // constructor
164 ~band_matrix(){}; // destructor
165 void
166 resize(int dim, int n_u, int n_l); // init with dim,n_u,n_l
167 int
168 dim() const; // matrix dimension
169 int
170 num_upper() const
171 {
172 return (int)d_upper.size() - 1;
173 }
174 int
175 num_lower() const
176 {
177 return (int)d_lower.size() - 1;
178 }
179 // access operator
180 double &
181 operator()(int i, int j); // write
182 double
183 operator()(int i, int j) const; // read
184 // we can store an additional diagonal (in d_lower)
185 double &
186 saved_diag(int i);
187 double
188 saved_diag(int i) const;
189 void
190 lu_decompose();
191 std::vector<double>
192 r_solve(const std::vector<double> &b) const;
193 std::vector<double>
194 l_solve(const std::vector<double> &b) const;
195 std::vector<double>
196 lu_solve(const std::vector<double> &b, bool is_lu_decomposed = false);
197 };
198 } // namespace splineInternal
199 } // namespace utils
200} // namespace dftefe
201
202#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:438
double get_x_max() const
Definition: Spline.h:142
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:117
double get_x_min() const
Definition: Spline.h:136
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:77
std::vector< double > get_y() const
Definition: Spline.h:131
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:368
bd_type d_right
Definition: Spline.h:71
std::vector< double > coefficients(double x) const
Definition: Spline.cpp:399
std::vector< double > d_d
Definition: Spline.h:68
bool d_made_monotonic
Definition: Spline.h:73
std::string info() const
Definition: Spline.cpp:501
size_t find_closest(double x) const
Definition: Spline.cpp:359
Spline()
Definition: Spline.cpp:44
std::vector< double > get_x() const
Definition: Spline.h:126
bool make_monotonic()
Definition: Spline.cpp:298
void set_coeffs_from_b()
Definition: Spline.cpp:91
std::vector< double > d_c
Definition: Spline.h:68
std::vector< double > d_x
Definition: Spline.h:64
std::vector< std::vector< double > > d_lower
Definition: Spline.h:160
int num_lower() const
Definition: Spline.h:175
band_matrix()
Definition: Spline.h:162
std::vector< std::vector< double > > d_upper
Definition: Spline.h:159
~band_matrix()
Definition: Spline.h:164
double & saved_diag(int i)
Definition: Spline.cpp:589
std::vector< double > r_solve(const std::vector< double > &b) const
Definition: Spline.cpp:658
int num_upper() const
Definition: Spline.h:170
std::vector< double > l_solve(const std::vector< double > &b) const
Definition: Spline.cpp:640
int dim() const
Definition: Spline.cpp:542
void resize(int dim, int n_u, int n_l)
Definition: Spline.cpp:525
double & operator()(int i, int j)
Definition: Spline.cpp:558
std::vector< double > lu_solve(const std::vector< double > &b, bool is_lu_decomposed=false)
Definition: Spline.cpp:676
void lu_decompose()
Definition: Spline.cpp:597
dealii includes
Definition: AtomFieldDataSpherical.cpp:31