ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
gsl::spline2d Class Reference

Higher level interface for interpolation. More...

#include <spline2d.hpp>

Collaboration diagram for gsl::spline2d:
Collaboration graph

Public Types

typedef gsl_interp2d_type type
 Convenience typedef. More...
 

Public Member Functions

 spline2d ()
 The default constructor is only really useful for assigning to. More...
 
 spline2d (type const *T, size_t const xsize, size_t const ysize)
 The default constructor creates a new spline2d with n elements. More...
 
 spline2d (gsl_spline2d *v)
 Could construct from a gsl_spline2d. More...
 
 spline2d (spline2d const &v)
 The copy constructor. More...
 
spline2doperator= (spline2d const &v)
 The assignment operator. More...
 
 ~spline2d ()
 The destructor only deletes the pointers if count reaches zero. More...
 
 spline2d (spline2d &&v)
 Move constructor. More...
 
spline2doperator= (spline2d &&v)
 Move operator. More...
 
bool operator== (spline2d const &v) const
 Two spline2d are identically equal if their elements are identical. More...
 
bool operator!= (spline2d const &v) const
 Two spline2d are different if their elements are not identical. More...
 
bool operator< (spline2d const &v) const
 A container needs to define an ordering for sorting. More...
 
bool operator> (spline2d const &v) const
 A container needs to define an ordering for sorting. More...
 
bool operator<= (spline2d const &v) const
 A container needs to define an ordering for sorting. More...
 
bool operator>= (spline2d const &v) const
 A container needs to define an ordering for sorting. More...
 
bool empty () const
 Find if the spline2d is empty. More...
 
void swap (spline2d &v)
 Swap two spline2d objects. More...
 
gsl_spline2d * get () const
 Get the gsl_spline2d. More...
 
bool unique () const
 Find if this is the only object sharing the gsl_spline2d. More...
 
size_t use_count () const
 Find how many spline2d objects share this pointer. More...
 
 operator bool () const
 Allow conversion to bool. More...
 
int init (double const xa[], double const ya[], double const za[], size_t xsize, size_t ysize)
 C++ version of gsl_spline2d_init(). More...
 
template<typename XA , typename YA , typename ZA >
int init (XA const &xa, YA const &ya, ZA const &za)
 C++ version of gsl_spline2d_init(). More...
 
double eval (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval(). More...
 
int eval_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_e(). More...
 
double eval_extrap (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval_extrap(). More...
 
int eval_extrap_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_extrap_e(). More...
 
double eval_deriv_x (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval_deriv_x(). More...
 
int eval_deriv_x_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_deriv_x_e(). More...
 
double eval_deriv_y (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval_deriv_y(). More...
 
int eval_deriv_y_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_deriv_y_e(). More...
 
double eval_deriv_xx (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval_deriv_xx(). More...
 
int eval_deriv_xx_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_deriv_xx_e(). More...
 
double eval_deriv_yy (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval_deriv_yy(). More...
 
int eval_deriv_yy_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_deriv_yy_e(). More...
 
double eval_deriv_xy (double const x, double const y, interp::accel &xa, interp::accel &ya) const
 C++ version of gsl_spline2d_eval_deriv_xy(). More...
 
int eval_deriv_xy_e (double const x, double const y, interp::accel &xa, interp::accel &ya, double &z) const
 C++ version of gsl_spline2d_eval_deriv_xy_e(). More...
 
size_t min_size () const
 C++ version of gsl_spline2d_min_size(). More...
 
char const * name () const
 C++ version of gsl_spline2d_name(). More...
 
int set (double zarr[], size_t const i, size_t const j, double const z) const
 C++ version of gsl_spline2d_set(). More...
 
template<typename ZARR >
int set (ZARR &zarr, size_t const i, size_t const j, double const z) const
 C++ version of gsl_spline2d_set(). More...
 
double get (double const zarr[], size_t const i, size_t const j) const
 C++ version of gsl_spline2d_get(). More...
 
template<typename ZARR >
double get (ZARR const &zarr, size_t const i, size_t const j) const
 C++ version of gsl_spline2d_get(). More...
 

Private Attributes

gsl_spline2d * ccgsl_pointer
 The shared pointer. More...
 
size_t * count
 The shared reference count. More...
 

Detailed Description

Higher level interface for interpolation.

Definition at line 32 of file spline2d.hpp.

Member Typedef Documentation

◆ type

typedef gsl_interp2d_type gsl::spline2d::type

Convenience typedef.

Definition at line 37 of file spline2d.hpp.

Constructor & Destructor Documentation

◆ spline2d() [1/5]

gsl::spline2d::spline2d ( )
inline

The default constructor is only really useful for assigning to.

Definition at line 41 of file spline2d.hpp.

References ccgsl_pointer, and count.

Referenced by operator=().

◆ spline2d() [2/5]

gsl::spline2d::spline2d ( type const *  T,
size_t const  xsize,
size_t const  ysize 
)
inlineexplicit

The default constructor creates a new spline2d with n elements.

Parameters
TThe interpolation type: use static functions from interp to choose a type
xsizeNumber of grid points in first direction
ysizeNumber of grid points in second direction

Definition at line 54 of file spline2d.hpp.

References ccgsl_pointer, and count.

◆ spline2d() [3/5]

gsl::spline2d::spline2d ( gsl_spline2d *  v)
inlineexplicit

Could construct from a gsl_spline2d.

This is not usually a good idea. In this case we should not use gsl_spline2d_free() to deallocate the memory.

Parameters
vThe spline2d

Definition at line 70 of file spline2d.hpp.

References ccgsl_pointer, and count.

◆ spline2d() [4/5]

gsl::spline2d::spline2d ( spline2d const &  v)
inline

The copy constructor.

This creates a new reference to the workspace.

Parameters
vThe spline2d to copy.

Definition at line 81 of file spline2d.hpp.

References ccgsl_pointer, and count.

◆ ~spline2d()

gsl::spline2d::~spline2d ( )
inline

The destructor only deletes the pointers if count reaches zero.

Definition at line 100 of file spline2d.hpp.

References ccgsl_pointer, and count.

◆ spline2d() [5/5]

gsl::spline2d::spline2d ( spline2d &&  v)
inline

Move constructor.

Parameters
vThe spline2d to move.

Definition at line 112 of file spline2d.hpp.

References count.

Member Function Documentation

◆ empty()

bool gsl::spline2d::empty ( ) const
inline

Find if the spline2d is empty.

Returns
true if has size zero; otherwise false

Definition at line 189 of file spline2d.hpp.

References ccgsl_pointer.

◆ eval()

double gsl::spline2d::eval ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
interpolated z value

Definition at line 268 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_x()

double gsl::spline2d::eval_deriv_x ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval_deriv_x().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
Interpolated derivative of z with respect to x

Definition at line 313 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_x_e()

int gsl::spline2d::eval_deriv_x_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_deriv_x_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zInterpolated derivative of z with respect to x
Returns
Error code on failure

Definition at line 325 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_xx()

double gsl::spline2d::eval_deriv_xx ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval_deriv_xx().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
Second derivative with respect to x

Definition at line 359 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_xx_e()

int gsl::spline2d::eval_deriv_xx_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_deriv_xx_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zSecond derivative with respect to x
Returns
Error code on failure

Definition at line 371 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_xy()

double gsl::spline2d::eval_deriv_xy ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval_deriv_xy().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
Second derivative with respect to x and y

Definition at line 405 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_xy_e()

int gsl::spline2d::eval_deriv_xy_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_deriv_xy_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zSecond derivative with respect to x and y
Returns
Error code on failure

Definition at line 417 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_y()

double gsl::spline2d::eval_deriv_y ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval_deriv_y().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
Interpolated derivative of z with respect to y

Definition at line 336 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_y_e()

int gsl::spline2d::eval_deriv_y_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_deriv_y_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zInterpolated derivative of z with respect to y
Returns
Error code on failure

Definition at line 348 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_yy()

double gsl::spline2d::eval_deriv_yy ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval_deriv_yy().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
Second derivative with respect to y

Definition at line 382 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_deriv_yy_e()

int gsl::spline2d::eval_deriv_yy_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_deriv_yy_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zSecond derivative with respect to y
Returns
Error code on failure

Definition at line 394 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_e()

int gsl::spline2d::eval_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zInterpolated or extrapolated z value
Returns
Error code on failure

Definition at line 279 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_extrap()

double gsl::spline2d::eval_extrap ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya 
) const
inline

C++ version of gsl_spline2d_eval_extrap().

Parameters
xx value
yy value
xax accelerator
yay accelerator
Returns
Interpolated or extrapolated z value

Definition at line 290 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ eval_extrap_e()

int gsl::spline2d::eval_extrap_e ( double const  x,
double const  y,
interp::accel xa,
interp::accel ya,
double &  z 
) const
inline

C++ version of gsl_spline2d_eval_extrap_e().

Parameters
xx value
yy value
xax accelerator
yay accelerator
zInterpolated or extrapolated z value [return]
Returns
error code on failure

Definition at line 302 of file spline2d.hpp.

References gsl::interp::accel::get(), and get().

◆ get() [1/3]

gsl_spline2d * gsl::spline2d::get ( ) const
inline

◆ get() [2/3]

double gsl::spline2d::get ( double const  zarr[],
size_t const  i,
size_t const  j 
) const
inline

C++ version of gsl_spline2d_get().

Parameters
zarr[]An array of size xsize × ysize
ix index in 0,…,xsize−1
jy index in 0,…,ysize−1
Returns
The value in zarr indexed by i and j

Definition at line 460 of file spline2d.hpp.

References get().

◆ get() [3/3]

template<typename ZARR >
double gsl::spline2d::get ( ZARR const &  zarr,
size_t const  i,
size_t const  j 
) const
inline

C++ version of gsl_spline2d_get().

This version handles std::vector and gsl::vector

Parameters
zarrAn array of size xsize × ysize
ix index in 0,…,xsize−1
jy index in 0,…,ysize−1
Returns
The value in zarr indexed by i and j

Definition at line 471 of file spline2d.hpp.

References get().

◆ init() [1/2]

int gsl::spline2d::init ( double const  xa[],
double const  ya[],
double const  za[],
size_t  xsize,
size_t  ysize 
)
inline

C++ version of gsl_spline2d_init().

Parameters
xa[]The array of x values
ya[]The array of y values
za[]The array of z values
xsizeThe size of xa
ysizeThe size of ya
Returns
Error code on failure

Definition at line 245 of file spline2d.hpp.

References get().

◆ init() [2/2]

template<typename XA , typename YA , typename ZA >
int gsl::spline2d::init ( XA const &  xa,
YA const &  ya,
ZA const &  za 
)
inline

C++ version of gsl_spline2d_init().

this version handles std::vector and gsl::vector.

Parameters
xaThe array of x values
yaThe array of y values
zaThe array of z values
Returns
Error code on failure

Definition at line 257 of file spline2d.hpp.

References get().

◆ min_size()

size_t gsl::spline2d::min_size ( ) const
inline

C++ version of gsl_spline2d_min_size().

Returns
minimum number of points required by this

Definition at line 424 of file spline2d.hpp.

References get().

◆ name()

char const * gsl::spline2d::name ( ) const
inline

C++ version of gsl_spline2d_name().

Returns
The name of the interpolation type

Definition at line 430 of file spline2d.hpp.

References get().

◆ operator bool()

gsl::spline2d::operator bool ( ) const
inlineexplicit

Allow conversion to bool.

Returns
true or false according as this contains a pointer to a gsl_spline2d.

Definition at line 235 of file spline2d.hpp.

References ccgsl_pointer.

◆ operator!=()

bool gsl::spline2d::operator!= ( spline2d const &  v) const
inline

Two spline2d are different if their elements are not identical.

Parameters
vThe spline2d to be compared with this
Returns
false or true according as this and v have identical elements or not

Definition at line 142 of file spline2d.hpp.

References operator==().

◆ operator<()

bool gsl::spline2d::operator< ( spline2d const &  v) const
inline

A container needs to define an ordering for sorting.

This uses standard lexicographical ordering and so is not useful, for example, for checking, that a spline2d is nonnegative.

Parameters
vThe spline2d to be compared with this
Returns
false or true according as this is less than v lexicographically

Definition at line 154 of file spline2d.hpp.

References ccgsl_pointer.

◆ operator<=()

bool gsl::spline2d::operator<= ( spline2d const &  v) const
inline

A container needs to define an ordering for sorting.

This uses standard lexicographical ordering and so is not useful, for example, for checking, that a spline2d is nonnegative.

Parameters
vThe spline2d to be compared with this
Returns
false or true according as this is less than or equal to v lexicographically

Definition at line 174 of file spline2d.hpp.

References ccgsl_pointer.

◆ operator=() [1/2]

spline2d & gsl::spline2d::operator= ( spline2d &&  v)
inline

Move operator.

Parameters
vThe spline2d to move.
Returns
A reference to this.

Definition at line 121 of file spline2d.hpp.

References spline2d().

◆ operator=() [2/2]

spline2d & gsl::spline2d::operator= ( spline2d const &  v)
inline

The assignment operator.

This copies elementwise.

Parameters
vThe spline2d to copy

Definition at line 88 of file spline2d.hpp.

References ccgsl_pointer, and count.

◆ operator==()

bool gsl::spline2d::operator== ( spline2d const &  v) const
inline

Two spline2d are identically equal if their elements are identical.

Parameters
vThe spline2d to be compared with this
Returns
true or false according as this and v have identical elements or not

Definition at line 134 of file spline2d.hpp.

References ccgsl_pointer.

Referenced by operator!=().

◆ operator>()

bool gsl::spline2d::operator> ( spline2d const &  v) const
inline

A container needs to define an ordering for sorting.

This uses standard lexicographical ordering and so is not useful, for example, for checking, that a spline2d is nonnegative.

Parameters
vThe spline2d to be compared with this
Returns
false or true according as this is greater than v lexicographically

Definition at line 164 of file spline2d.hpp.

References ccgsl_pointer.

◆ operator>=()

bool gsl::spline2d::operator>= ( spline2d const &  v) const
inline

A container needs to define an ordering for sorting.

This uses standard lexicographical ordering and so is not useful, for example, for checking, that a spline2d is nonnegative.

Parameters
vThe spline2d to be compared with this
Returns
false or true according as this is no less than v lexicographically

Definition at line 184 of file spline2d.hpp.

References ccgsl_pointer.

◆ set() [1/2]

int gsl::spline2d::set ( double  zarr[],
size_t const  i,
size_t const  j,
double const  z 
) const
inline

C++ version of gsl_spline2d_set().

Parameters
zarr[]An array of size xsize × ysize
ix index in 0,…,xsize−1
jy index in 0,…,ysize−1
zThe value to be set
Returns
Error code on failure

Definition at line 439 of file spline2d.hpp.

References get().

◆ set() [2/2]

template<typename ZARR >
int gsl::spline2d::set ( ZARR &  zarr,
size_t const  i,
size_t const  j,
double const  z 
) const
inline

C++ version of gsl_spline2d_set().

This version handles std::vector and gsl::vector

Parameters
zarrAn array of size xsize × ysize
ix index in 0,…,xsize−1
jy index in 0,…,ysize−1
zThe value to be set
Returns
Error code on failure

Definition at line 451 of file spline2d.hpp.

References get().

◆ swap()

void gsl::spline2d::swap ( spline2d v)
inline

Swap two spline2d objects.

This works even if the spline2d have different sizes because it swaps pointers.

Parameters
vThe spline2d to swap with this.

Definition at line 196 of file spline2d.hpp.

References ccgsl_pointer, and count.

◆ unique()

bool gsl::spline2d::unique ( ) const
inline

Find if this is the only object sharing the gsl_spline2d.

Returns
true or falses according as this is the only spline2d object sharing the gsl_spline2d.

Definition at line 221 of file spline2d.hpp.

References count.

◆ use_count()

size_t gsl::spline2d::use_count ( ) const
inline

Find how many spline2d objects share this pointer.

Returns
the number of spline2d objects that share this pointer.

Definition at line 226 of file spline2d.hpp.

References count.

Member Data Documentation

◆ ccgsl_pointer

gsl_spline2d* gsl::spline2d::ccgsl_pointer
private

The shared pointer.

Definition at line 204 of file spline2d.hpp.

Referenced by empty(), get(), operator bool(), operator<(), operator<=(), operator=(), operator==(), operator>(), operator>=(), spline2d(), swap(), and ~spline2d().

◆ count

size_t* gsl::spline2d::count
private

The shared reference count.

Definition at line 208 of file spline2d.hpp.

Referenced by operator=(), spline2d(), swap(), unique(), use_count(), and ~spline2d().


The documentation for this class was generated from the following file: