ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
spline2d.hpp
Go to the documentation of this file.
1/*
2 * $Id$
3 * Copyright (C) 2012, 2019, 2021 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_SPLINE2D_HPP
21#define CCGSL_SPLINE2D_HPP
22
23#include<new>
24#include<gsl/gsl_spline2d.h>
25#include"interp.hpp"
26#include"interp2d.hpp"
27
28namespace gsl {
32 class spline2d {
33 public:
37 typedef gsl_interp2d_type type;
42 ccgsl_pointer = 0;
43 count = 0; // initially nullptr will do
44 }
45 // Refines random access container
46 // Refines assignable
54 explicit spline2d( type const* T, size_t const xsize, size_t const ysize ){
55 ccgsl_pointer = gsl_spline2d_alloc( T, xsize, ysize );
56 // just plausibly we could allocate spline2d but not count
57 try { count = new size_t; } catch( std::bad_alloc& e ){
58 // try to tidy up before rethrowing
59 gsl_spline2d_free( ccgsl_pointer );
60 throw e;
61 }
62 *count = 1; // initially there is just one reference to ccgsl_pointer
63 }
70 explicit spline2d( gsl_spline2d* v ){
71 ccgsl_pointer = v;
72 // just plausibly we could fail to allocate count: no further action needed.
73 count = new size_t;
74 *count = 1; // initially there is just one reference to ccgsl_pointer
75 }
76 // copy constructor
82 count = v.count; if( count != 0 ) ++*count; }
83 // assignment operator
89 // first, possibly delete anything pointed to by this
90 if( count == 0 or --*count == 0 ){
91 if( ccgsl_pointer != 0 ) gsl_spline2d_free( ccgsl_pointer );
92 delete count;
93 } // Then copy
94 ccgsl_pointer = v.ccgsl_pointer; count = v.count; if( count != 0 ) ++*count; return *this;
95 }
96 // destructor
101 if( count == 0 or --*count == 0 ){
102 // could have allocated null pointer
103 if( ccgsl_pointer != 0 ) gsl_spline2d_free( ccgsl_pointer );
104 delete count;
105 }
106 }
107#ifdef __GXX_EXPERIMENTAL_CXX0X__
113 std::swap( count, v.count );
114 v.ccgsl_pointer = nullptr;
115 }
122 spline2d( std::move( v ) ).swap( *this );
123 return *this;
124 }
125#endif
126 // Refines equality comparable
127 // == operator
134 bool operator==( spline2d const& v ) const { return ccgsl_pointer == v.ccgsl_pointer; }
135 // != operator
142 bool operator!=( spline2d const& v ) const { return not operator==( v ); }
143 // Refines forward container
144 // Refines less than comparable
145 // operator<
154 bool operator<( spline2d const& v ) const { return ccgsl_pointer < v.ccgsl_pointer; }
155 // operator>
164 bool operator>( spline2d const& v ) const { return ccgsl_pointer > v.ccgsl_pointer; }
165 // operator<=
174 bool operator<=( spline2d const& v ) const { return ccgsl_pointer <= v.ccgsl_pointer; }
175 // operator>=
184 bool operator>=( spline2d const& v ) const { return ccgsl_pointer >= v.ccgsl_pointer; }
189 bool empty() const { return ccgsl_pointer == 0; }
190 // swap() --- should work even if sizes don't match
196 void swap( spline2d& v ){
197 std::swap( ccgsl_pointer, v.ccgsl_pointer );
198 std::swap( count, v.count );
199 }
200 private:
204 gsl_spline2d* ccgsl_pointer;
208 size_t* count;
209 public:
210 // shared reference functions
215 gsl_spline2d* get() const { return ccgsl_pointer; }
221 bool unique() const { return count != 0 and *count == 1; }
226 size_t use_count() const { return count == 0 ? 0 : *count; }
232#ifdef __GXX_EXPERIMENTAL_CXX0X__
233 explicit
234#endif
235 operator bool() const { return ccgsl_pointer != 0; }
245 int init( double const xa[], double const ya[], double const za[],
246 size_t xsize, size_t ysize ){
247 return gsl_spline2d_init( get(), xa, ya, za, xsize, ysize ); }
256 template<typename XA,typename YA, typename ZA>
257 int init( XA const& xa, YA const& ya, ZA const& za ){
258 return gsl_spline2d_init( get(), xa.data(), ya.data(), za.data(),
259 xa.size(), ya.size() ); }
268 double eval( double const x, double const y, interp::accel& xa, interp::accel& ya ) const {
269 return gsl_spline2d_eval( get(), x, y, xa.get(), ya.get() ); }
279 int eval_e( double const x, double const y, interp::accel& xa, interp::accel& ya,
280 double& z ) const {
281 return gsl_spline2d_eval_e( get(), x, y, xa.get(), ya.get(), &z ); }
290 double eval_extrap( double const x, double const y, interp::accel& xa,
291 interp::accel& ya ) const {
292 return gsl_spline2d_eval_extrap( get(), x, y, xa.get(), ya.get() ); }
302 int eval_extrap_e( double const x, double const y,
303 interp::accel& xa, interp::accel& ya, double& z ) const {
304 return gsl_spline2d_eval_extrap_e( get(), x, y, xa.get(), ya.get(), &z ); }
313 double eval_deriv_x( double const x, double const y, interp::accel& xa,
314 interp::accel& ya ) const {
315 return gsl_spline2d_eval_deriv_x( get(), x, y, xa.get(), ya.get() ); }
325 int eval_deriv_x_e( double const x, double const y, interp::accel& xa, interp::accel& ya,
326 double& z ) const {
327 return gsl_spline2d_eval_deriv_x_e( get(), x, y, xa.get(), ya.get(), &z ); }
336 double eval_deriv_y( double const x, double const y, interp::accel& xa,
337 interp::accel& ya ) const {
338 return gsl_spline2d_eval_deriv_y( get(), x, y, xa.get(), ya.get() ); }
348 int eval_deriv_y_e( double const x, double const y, interp::accel& xa, interp::accel& ya,
349 double& z ) const {
350 return gsl_spline2d_eval_deriv_y_e( get(), x, y, xa.get(), ya.get(), &z ); }
359 double eval_deriv_xx( double const x, double const y, interp::accel& xa,
360 interp::accel& ya ) const {
361 return gsl_spline2d_eval_deriv_xx( get(), x, y, xa.get(), ya.get() ); }
371 int eval_deriv_xx_e( double const x, double const y, interp::accel& xa, interp::accel& ya,
372 double& z ) const {
373 return gsl_spline2d_eval_deriv_xx_e( get(), x, y, xa.get(), ya.get(), &z ); }
382 double eval_deriv_yy( double const x, double const y, interp::accel& xa,
383 interp::accel& ya ) const {
384 return gsl_spline2d_eval_deriv_yy( get(), x, y, xa.get(), ya.get() ); }
394 int eval_deriv_yy_e( double const x, double const y, interp::accel& xa, interp::accel& ya,
395 double& z ) const {
396 return gsl_spline2d_eval_deriv_yy_e( get(), x, y, xa.get(), ya.get(), &z ); }
405 double eval_deriv_xy( double const x, double const y, interp::accel& xa,
406 interp::accel& ya ) const {
407 return gsl_spline2d_eval_deriv_xy( get(), x, y, xa.get(), ya.get() ); }
417 int eval_deriv_xy_e( double const x, double const y, interp::accel& xa, interp::accel& ya,
418 double& z ) const {
419 return gsl_spline2d_eval_deriv_xy_e( get(), x, y, xa.get(), ya.get(), &z ); }
424 size_t min_size() const { return gsl_spline2d_min_size( get() ); }
429 char const*
430 name() const { return gsl_spline2d_name( get() ); }
439 int set( double zarr[], size_t const i, size_t const j, double const z ) const {
440 return gsl_spline2d_set( get(), zarr, i, j, z ); }
450 template<typename ZARR>
451 int set( ZARR& zarr, size_t const i, size_t const j, double const z ) const {
452 return gsl_spline2d_set( get(), zarr.data(), i, j, z ); }
460 double get( double const zarr[], size_t const i, size_t const j ) const {
461 return gsl_spline2d_get( get(), zarr, i, j ); }
470 template<typename ZARR>
471 double get( ZARR const& zarr, size_t const i, size_t const j ) const {
472 return gsl_spline2d_get( get(), zarr.data(), i, j ); }
473 };
474}
475#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: spline2d.hpp:32
gsl_spline2d * ccgsl_pointer
The shared pointer.
Definition: spline2d.hpp:204
size_t min_size() const
C++ version of gsl_spline2d_min_size().
Definition: spline2d.hpp:424
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().
Definition: spline2d.hpp:336
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().
Definition: spline2d.hpp:417
int set(ZARR &zarr, size_t const i, size_t const j, double const z) const
C++ version of gsl_spline2d_set().
Definition: spline2d.hpp:451
spline2d & operator=(spline2d const &v)
The assignment operator.
Definition: spline2d.hpp:88
int set(double zarr[], size_t const i, size_t const j, double const z) const
C++ version of gsl_spline2d_set().
Definition: spline2d.hpp:439
int init(double const xa[], double const ya[], double const za[], size_t xsize, size_t ysize)
C++ version of gsl_spline2d_init().
Definition: spline2d.hpp:245
size_t * count
The shared reference count.
Definition: spline2d.hpp:208
spline2d(spline2d &&v)
Move constructor.
Definition: spline2d.hpp:112
spline2d(spline2d const &v)
The copy constructor.
Definition: spline2d.hpp:81
bool operator!=(spline2d const &v) const
Two spline2d are different if their elements are not identical.
Definition: spline2d.hpp:142
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().
Definition: spline2d.hpp:279
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().
Definition: spline2d.hpp:348
gsl_interp2d_type type
Convenience typedef.
Definition: spline2d.hpp:37
double eval_extrap(double const x, double const y, interp::accel &xa, interp::accel &ya) const
C++ version of gsl_spline2d_eval_extrap().
Definition: spline2d.hpp:290
spline2d & operator=(spline2d &&v)
Move operator.
Definition: spline2d.hpp:121
size_t use_count() const
Find how many spline2d objects share this pointer.
Definition: spline2d.hpp:226
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().
Definition: spline2d.hpp:371
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().
Definition: spline2d.hpp:382
bool operator<(spline2d const &v) const
A container needs to define an ordering for sorting.
Definition: spline2d.hpp:154
spline2d(gsl_spline2d *v)
Could construct from a gsl_spline2d.
Definition: spline2d.hpp:70
double get(ZARR const &zarr, size_t const i, size_t const j) const
C++ version of gsl_spline2d_get().
Definition: spline2d.hpp:471
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().
Definition: spline2d.hpp:359
void swap(spline2d &v)
Swap two spline2d objects.
Definition: spline2d.hpp:196
gsl_spline2d * get() const
Get the gsl_spline2d.
Definition: spline2d.hpp:215
int init(XA const &xa, YA const &ya, ZA const &za)
C++ version of gsl_spline2d_init().
Definition: spline2d.hpp:257
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().
Definition: spline2d.hpp:405
bool operator==(spline2d const &v) const
Two spline2d are identically equal if their elements are identical.
Definition: spline2d.hpp:134
bool unique() const
Find if this is the only object sharing the gsl_spline2d.
Definition: spline2d.hpp:221
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().
Definition: spline2d.hpp:302
bool empty() const
Find if the spline2d is empty.
Definition: spline2d.hpp:189
char const * name() const
C++ version of gsl_spline2d_name().
Definition: spline2d.hpp:430
bool operator>(spline2d const &v) const
A container needs to define an ordering for sorting.
Definition: spline2d.hpp:164
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().
Definition: spline2d.hpp:394
spline2d()
The default constructor is only really useful for assigning to.
Definition: spline2d.hpp:41
bool operator>=(spline2d const &v) const
A container needs to define an ordering for sorting.
Definition: spline2d.hpp:184
spline2d(type const *T, size_t const xsize, size_t const ysize)
The default constructor creates a new spline2d with n elements.
Definition: spline2d.hpp:54
double get(double const zarr[], size_t const i, size_t const j) const
C++ version of gsl_spline2d_get().
Definition: spline2d.hpp:460
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().
Definition: spline2d.hpp:313
~spline2d()
The destructor only deletes the pointers if count reaches zero.
Definition: spline2d.hpp:100
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().
Definition: spline2d.hpp:325
bool operator<=(spline2d const &v) const
A container needs to define an ordering for sorting.
Definition: spline2d.hpp:174
double eval(double const x, double const y, interp::accel &xa, interp::accel &ya) const
C++ version of gsl_spline2d_eval().
Definition: spline2d.hpp:268
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34