20#ifndef CCGSL_SPLINE_HPP
21#define CCGSL_SPLINE_HPP
24#include<gsl/gsl_spline.h>
36 typedef gsl_interp_type
type;
55 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
105#ifdef __GXX_EXPERIMENTAL_CXX0X__
111 std::swap(
count, v.count );
112 v.ccgsl_pointer =
nullptr;
120 spline( std::move( v ) ).swap( *
this );
230#ifdef __GXX_EXPERIMENTAL_CXX0X__
242 int init(
double const xa[],
double const ya[],
size_t size ){
243 return gsl_spline_init(
get(), xa, ya,
size ); }
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() ); }
260 char const*
name()
const {
return gsl_spline_name(
get() ); }
266 unsigned int min_size()
const {
return gsl_spline_min_size(
get() ); }
276 return gsl_spline_eval_e(
get(), x,
a.get(), y ); }
285 return gsl_spline_eval_e(
get(), x,
a.get(), &y ); }
294 return gsl_spline_eval(
get(), x,
a.get() ); }
304 return gsl_spline_eval_deriv_e(
get(), x,
a.get(), y ); }
313 return gsl_spline_eval_deriv_e(
get(), x,
a.get(), &y ); }
322 return gsl_spline_eval_deriv(
get(), x,
a.get() ); }
332 return gsl_spline_eval_deriv2_e(
get(), x,
a.get(), y ); }
341 return gsl_spline_eval_deriv2_e(
get(), x,
a.get(), &y ); }
350 return gsl_spline_eval_deriv2(
get(), x,
a.get() ); }
361 return gsl_spline_eval_integ_e(
get(),
a,
b, acc.
get(), y ); }
371 return gsl_spline_eval_integ_e(
get(),
a,
b, acc.
get(), &y ); }
381 return gsl_spline_eval_integ(
get(),
a,
b, acc.
get() ); }
Workspace for acceleration.
gsl_interp_accel * get() const
Get the gsl_interp_accel.
Higher level interface for interpolation.
bool operator<=(spline const &v) const
A container needs to define an ordering for sorting.
spline & operator=(spline &&v)
Move operator.
spline(gsl_spline *v)
Could construct from a gsl_spline.
bool operator>=(spline const &v) const
A container needs to define an ordering for sorting.
double eval_deriv2(double x, interp::accel &a) const
C++ version of gsl_spline_eval_deriv2().
bool operator==(spline const &v) const
Two spline are identically equal if their elements are identical.
void swap(spline &v)
Swap two spline objects.
int init(double const xa[], double const ya[], size_t size)
C++ version of gsl_spline_init().
int eval_e(double x, interp::accel &a, double &y) const
C++ version of gsl_spline_eval_e().
double eval_deriv(double x, interp::accel &a) const
C++ version of gsl_spline_eval_deriv().
size_t use_count() const
Find how many spline objects share this pointer.
size_t * count
The shared reference count.
int eval_deriv2_e(double x, interp::accel &a, double *y) const
C++ version of gsl_spline_eval_deriv2_e().
bool operator!=(spline const &v) const
Two spline are different if their elements are not identical.
int eval_integ_e(double a, double b, interp::accel &acc, double &y) const
C++ version of gsl_spline_eval_integ_e().
gsl_interp_type type
Convenience typedef.
unsigned int min_size() const
C++ version of gsl_spline_min_size().
spline(spline const &v)
The copy constructor.
bool empty() const
Find if the spline is empty.
int init(XA const &xa, YA const &ya)
C++ version of gsl_spline_init().
int eval_integ_e(double a, double b, interp::accel &acc, double *y) const
C++ version of gsl_spline_eval_integ_e().
double eval_integ(double a, double b, interp::accel &acc) const
C++ version of gsl_spline_eval_integ().
bool operator>(spline const &v) const
A container needs to define an ordering for sorting.
char const * name() const
C++ version of gsl_spline_name().
int eval_deriv_e(double x, interp::accel &a, double *y) const
C++ version of gsl_spline_eval_deriv_e().
bool operator<(spline const &v) const
A container needs to define an ordering for sorting.
int eval_deriv2_e(double x, interp::accel &a, double &y) const
C++ version of gsl_spline_eval_deriv2_e().
~spline()
The destructor only deletes the pointers if count reaches zero.
int eval_e(double x, interp::accel &a, double *y) const
C++ version of gsl_spline_eval_e().
spline()
The default constructor is only really useful for assigning to.
gsl_spline * ccgsl_pointer
The shared pointer.
bool unique() const
Find if this is the only object sharing the gsl_spline.
spline(type const *T, size_t const n)
The default constructor creates a new spline with n elements.
spline(spline &&v)
Move constructor.
double eval(double x, interp::accel &a) const
C++ version of gsl_spline_eval().
spline & operator=(spline const &v)
The assignment operator.
int eval_deriv_e(double x, interp::accel &a, double &y) const
C++ version of gsl_spline_eval_deriv_e().
gsl_spline * get() const
Get the gsl_spline.
size_t size(series const &cs)
C++ version of gsl_cheb_size().
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
double a(int order, double qq)
C++ version of gsl_sf_mathieu_a().
The gsl package creates an interface to the GNU Scientific Library for C++.