ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
spline.hpp
Go to the documentation of this file.
1/*
2 * $Id: spline.hpp 293 2012-12-17 20:27:36Z jdl3 $
3 * Copyright (C) 2012 John D Lamb
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or (at
8 * your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18 */
19
20#ifndef CCGSL_SPLINE_HPP
21#define CCGSL_SPLINE_HPP
22
23#include<new>
24#include<gsl/gsl_spline.h>
25#include"interp.hpp"
26
27namespace gsl {
31 class spline {
32 public:
36 typedef gsl_interp_type type;
41 ccgsl_pointer = 0;
42 count = 0; // initially nullptr will do
43 }
44 // Refines random access container
45 // Refines assignable
52 explicit spline( type const* T, size_t const n ){
53 ccgsl_pointer = gsl_spline_alloc( T, n );
54 // just plausibly we could allocate spline but not count
55 try { count = new size_t; } catch( std::bad_alloc& e ){
56 // try to tidy up before rethrowing
57 gsl_spline_free( ccgsl_pointer );
58 throw e;
59 }
60 *count = 1; // initially there is just one reference to ccgsl_pointer
61 }
68 explicit spline( gsl_spline* v ){
69 ccgsl_pointer = v;
70 // just plausibly we could fail to allocate count: no further action needed.
71 count = new size_t;
72 *count = 1; // initially there is just one reference to ccgsl_pointer
73 }
74 // copy constructor
80 count = v.count; if( count != 0 ) ++*count; }
81 // assignment operator
86 spline& operator=( spline const& v ){
87 // first, possibly delete anything pointed to by this
88 if( count == 0 or --*count == 0 ){
89 if( ccgsl_pointer != 0 ) gsl_spline_free( ccgsl_pointer );
90 delete count;
91 } // Then copy
92 ccgsl_pointer = v.ccgsl_pointer; count = v.count; if( count != 0 ) ++*count; return *this;
93 }
94 // destructor
99 if( count == 0 or --*count == 0 ){
100 // could have allocated null pointer
101 if( ccgsl_pointer != 0 ) gsl_spline_free( ccgsl_pointer );
102 delete count;
103 }
104 }
105#ifdef __GXX_EXPERIMENTAL_CXX0X__
110 spline( spline&& v ) : ccgsl_pointer( v.ccgsl_pointer ), count( nullptr ){
111 std::swap( count, v.count );
112 v.ccgsl_pointer = nullptr;
113 }
120 spline( std::move( v ) ).swap( *this );
121 return *this;
122 }
123#endif
124 // Refines equality comparable
125 // == operator
132 bool operator==( spline const& v ) const { return ccgsl_pointer == v.ccgsl_pointer; }
133 // != operator
140 bool operator!=( spline const& v ) const { return not operator==( v ); }
141 // Refines forward container
142 // Refines less than comparable
143 // operator<
152 bool operator<( spline const& v ) const { return ccgsl_pointer < v.ccgsl_pointer; }
153 // operator>
162 bool operator>( spline const& v ) const { return ccgsl_pointer > v.ccgsl_pointer; }
163 // operator<=
172 bool operator<=( spline const& v ) const { return ccgsl_pointer <= v.ccgsl_pointer; }
173 // operator>=
182 bool operator>=( spline const& v ) const { return ccgsl_pointer >= v.ccgsl_pointer; }
187 bool empty() const { return ccgsl_pointer == 0; }
188 // swap() --- should work even if sizes don't match
194 void swap( spline& v ){
195 std::swap( ccgsl_pointer, v.ccgsl_pointer );
196 std::swap( count, v.count );
197 }
198 private:
202 gsl_spline* ccgsl_pointer;
206 size_t* count;
207 public:
208 // shared reference functions
213 gsl_spline* get() const { return ccgsl_pointer; }
219 bool unique() const { return count != 0 and *count == 1; }
224 size_t use_count() const { return count == 0 ? 0 : *count; }
230#ifdef __GXX_EXPERIMENTAL_CXX0X__
231 explicit
232#endif
233 operator bool() const { return ccgsl_pointer != 0; }
234
242 int init( double const xa[], double const ya[], size_t size ){
243 return gsl_spline_init( get(), xa, ya, size ); }
244
252 template<typename XA,typename YA>
253 int init( XA const& xa, YA const& ya ){
254 return gsl_spline_init( get(), xa.data(), ya.data(), xa.size() ); }
255
260 char const* name() const { return gsl_spline_name( get() ); }
261
266 unsigned int min_size() const { return gsl_spline_min_size( get() ); }
267
275 int eval_e( double x, interp::accel& a, double* y ) const {
276 return gsl_spline_eval_e( get(), x, a.get(), y ); }
284 int eval_e( double x, interp::accel& a, double& y ) const {
285 return gsl_spline_eval_e( get(), x, a.get(), &y ); }
286
293 double eval( double x, interp::accel& a ) const {
294 return gsl_spline_eval( get(), x, a.get() ); }
295
303 int eval_deriv_e( double x, interp::accel& a, double* y ) const {
304 return gsl_spline_eval_deriv_e( get(), x, a.get(), y ); }
312 int eval_deriv_e( double x, interp::accel& a, double& y ) const {
313 return gsl_spline_eval_deriv_e( get(), x, a.get(), &y ); }
314
321 double eval_deriv( double x, interp::accel& a ) const {
322 return gsl_spline_eval_deriv( get(), x, a.get() ); }
323
331 int eval_deriv2_e( double x, interp::accel& a, double* y ) const {
332 return gsl_spline_eval_deriv2_e( get(), x, a.get(), y ); }
340 int eval_deriv2_e( double x, interp::accel& a, double& y ) const {
341 return gsl_spline_eval_deriv2_e( get(), x, a.get(), &y ); }
342
349 double eval_deriv2( double x, interp::accel& a ) const {
350 return gsl_spline_eval_deriv2( get(), x, a.get() ); }
351
360 int eval_integ_e( double a, double b, interp::accel& acc, double* y ) const {
361 return gsl_spline_eval_integ_e( get(), a, b, acc.get(), y ); }
370 int eval_integ_e( double a, double b, interp::accel& acc, double& y ) const {
371 return gsl_spline_eval_integ_e( get(), a, b, acc.get(), &y ); }
372
380 double eval_integ( double a, double b, interp::accel& acc ) const {
381 return gsl_spline_eval_integ( get(), a, b, acc.get() ); }
382
383 };
384}
385#endif
Workspace for acceleration.
Definition: interp.hpp:237
gsl_interp_accel * get() const
Get the gsl_interp_accel.
Definition: interp.hpp:405
Higher level interface for interpolation.
Definition: spline.hpp:31
bool operator<=(spline const &v) const
A container needs to define an ordering for sorting.
Definition: spline.hpp:172
spline & operator=(spline &&v)
Move operator.
Definition: spline.hpp:119
spline(gsl_spline *v)
Could construct from a gsl_spline.
Definition: spline.hpp:68
bool operator>=(spline const &v) const
A container needs to define an ordering for sorting.
Definition: spline.hpp:182
double eval_deriv2(double x, interp::accel &a) const
C++ version of gsl_spline_eval_deriv2().
Definition: spline.hpp:349
bool operator==(spline const &v) const
Two spline are identically equal if their elements are identical.
Definition: spline.hpp:132
void swap(spline &v)
Swap two spline objects.
Definition: spline.hpp:194
int init(double const xa[], double const ya[], size_t size)
C++ version of gsl_spline_init().
Definition: spline.hpp:242
int eval_e(double x, interp::accel &a, double &y) const
C++ version of gsl_spline_eval_e().
Definition: spline.hpp:284
double eval_deriv(double x, interp::accel &a) const
C++ version of gsl_spline_eval_deriv().
Definition: spline.hpp:321
size_t use_count() const
Find how many spline objects share this pointer.
Definition: spline.hpp:224
size_t * count
The shared reference count.
Definition: spline.hpp:206
int eval_deriv2_e(double x, interp::accel &a, double *y) const
C++ version of gsl_spline_eval_deriv2_e().
Definition: spline.hpp:331
bool operator!=(spline const &v) const
Two spline are different if their elements are not identical.
Definition: spline.hpp:140
int eval_integ_e(double a, double b, interp::accel &acc, double &y) const
C++ version of gsl_spline_eval_integ_e().
Definition: spline.hpp:370
gsl_interp_type type
Convenience typedef.
Definition: spline.hpp:36
unsigned int min_size() const
C++ version of gsl_spline_min_size().
Definition: spline.hpp:266
spline(spline const &v)
The copy constructor.
Definition: spline.hpp:79
bool empty() const
Find if the spline is empty.
Definition: spline.hpp:187
int init(XA const &xa, YA const &ya)
C++ version of gsl_spline_init().
Definition: spline.hpp:253
int eval_integ_e(double a, double b, interp::accel &acc, double *y) const
C++ version of gsl_spline_eval_integ_e().
Definition: spline.hpp:360
double eval_integ(double a, double b, interp::accel &acc) const
C++ version of gsl_spline_eval_integ().
Definition: spline.hpp:380
bool operator>(spline const &v) const
A container needs to define an ordering for sorting.
Definition: spline.hpp:162
char const * name() const
C++ version of gsl_spline_name().
Definition: spline.hpp:260
int eval_deriv_e(double x, interp::accel &a, double *y) const
C++ version of gsl_spline_eval_deriv_e().
Definition: spline.hpp:303
bool operator<(spline const &v) const
A container needs to define an ordering for sorting.
Definition: spline.hpp:152
int eval_deriv2_e(double x, interp::accel &a, double &y) const
C++ version of gsl_spline_eval_deriv2_e().
Definition: spline.hpp:340
~spline()
The destructor only deletes the pointers if count reaches zero.
Definition: spline.hpp:98
int eval_e(double x, interp::accel &a, double *y) const
C++ version of gsl_spline_eval_e().
Definition: spline.hpp:275
spline()
The default constructor is only really useful for assigning to.
Definition: spline.hpp:40
gsl_spline * ccgsl_pointer
The shared pointer.
Definition: spline.hpp:202
bool unique() const
Find if this is the only object sharing the gsl_spline.
Definition: spline.hpp:219
spline(type const *T, size_t const n)
The default constructor creates a new spline with n elements.
Definition: spline.hpp:52
spline(spline &&v)
Move constructor.
Definition: spline.hpp:110
double eval(double x, interp::accel &a) const
C++ version of gsl_spline_eval().
Definition: spline.hpp:293
spline & operator=(spline const &v)
The assignment operator.
Definition: spline.hpp:86
int eval_deriv_e(double x, interp::accel &a, double &y) const
C++ version of gsl_spline_eval_deriv_e().
Definition: spline.hpp:312
gsl_spline * get() const
Get the gsl_spline.
Definition: spline.hpp:213
size_t size(series const &cs)
C++ version of gsl_cheb_size().
Definition: chebyshev.hpp:287
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
Definition: rstat.hpp:299
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
Definition: sf_mathieu.hpp:298
double a(int order, double qq)
C++ version of gsl_sf_mathieu_a().
Definition: sf_mathieu.hpp:272
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34