BubbleProfiler  0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
BubbleProfiler::GSL_interpolator Class Reference

Provides routines for carrying out 1D interpolation. More...

#include <gsl_interpolator.hpp>

Classes

struct  GSL_accel_deleter
 
struct  GSL_spline_deleter
 

Public Types

enum  Interpolation_type {
  Interpolation_type::GSL_linear, Interpolation_type::GSL_polynomial, Interpolation_type::GSL_cspline, Interpolation_type::GSL_cspline_periodic,
  Interpolation_type::GSL_akima, Interpolation_type::GSL_akima_periodic
}
 

Public Member Functions

 GSL_interpolator ()
 
 ~GSL_interpolator ()=default
 
 GSL_interpolator (const GSL_interpolator &)
 
GSL_interpolatoroperator= (const GSL_interpolator &)
 
 GSL_interpolator (GSL_interpolator &&)=default
 
GSL_interpolatoroperator= (GSL_interpolator &&)=default
 
void set_interpolation_type (Interpolation_type t)
 Sets the method used to perform interpolation. More...
 
void construct_interpolant (const Eigen::VectorXd &, const Eigen::VectorXd &)
 Computes interpolant using given x and y values. More...
 
void construct_interpolant (const Eigen::VectorXd &, const Eigen::VectorXd &, int)
 Computes interpolant using requested number of points. More...
 
double evaluate_at (double) const
 Evaluates computed interpolant at a point. More...
 
double evaluate_deriv_at (double) const
 Evaluates derivative of computed interpolant at a point. More...
 
double evaluate_second_deriv_at (double) const
 Evaluates second derivative of computed interpolant at a point. More...
 

Private Member Functions

const gsl_interp_type * to_gsl_interp_type (Interpolation_type) const
 
void reset (Interpolation_type, int)
 
void check_interpolant_computed () const
 
void check_within_bounds (double) const
 

Private Attributes

Interpolation_type interp_type {Interpolation_type::GSL_cspline}
 
std::unique_ptr< gsl_interp_accel, GSL_accel_deleteraccelerator {nullptr}
 
std::unique_ptr< gsl_spline, GSL_spline_deleterinterpolant {nullptr}
 

Detailed Description

Provides routines for carrying out 1D interpolation.

Definition at line 64 of file gsl_interpolator.hpp.

Member Enumeration Documentation

Enumerator
GSL_linear 
GSL_polynomial 
GSL_cspline 
GSL_cspline_periodic 
GSL_akima 
GSL_akima_periodic 

Definition at line 66 of file gsl_interpolator.hpp.

Constructor & Destructor Documentation

BubbleProfiler::GSL_interpolator::GSL_interpolator ( )

Definition at line 176 of file gsl_interpolator.cpp.

References accelerator.

BubbleProfiler::GSL_interpolator::~GSL_interpolator ( )
default
BubbleProfiler::GSL_interpolator::GSL_interpolator ( const GSL_interpolator other)

Definition at line 185 of file gsl_interpolator.cpp.

References accelerator, interp_type, and interpolant.

BubbleProfiler::GSL_interpolator::GSL_interpolator ( GSL_interpolator &&  )
default

Member Function Documentation

void BubbleProfiler::GSL_interpolator::check_interpolant_computed ( ) const
private

Definition at line 414 of file gsl_interpolator.cpp.

References interpolant.

Referenced by evaluate_at(), evaluate_deriv_at(), and evaluate_second_deriv_at().

void BubbleProfiler::GSL_interpolator::check_within_bounds ( double  x) const
private

Definition at line 422 of file gsl_interpolator.cpp.

References interpolant.

Referenced by evaluate_at(), evaluate_deriv_at(), and evaluate_second_deriv_at().

void BubbleProfiler::GSL_interpolator::construct_interpolant ( const Eigen::VectorXd &  x,
const Eigen::VectorXd &  y 
)

Computes interpolant using given x and y values.

Parameters
[in]xthe x values to use
[in]ythe function values to use

The interpolant is computed using all of the given (x, y) values. The minimum required number of points depends on the choice of interpolation method. The vectors x and y must have the same lengths. The resulting interpolant can subsequently be evaluated by calls to evaluate_at , evaluate_deriv_at , and evaluate_second_deriv_at .

Definition at line 233 of file gsl_interpolator.cpp.

Referenced by set_interpolation_type().

void BubbleProfiler::GSL_interpolator::construct_interpolant ( const Eigen::VectorXd &  x,
const Eigen::VectorXd &  y,
int  n 
)

Computes interpolant using requested number of points.

Parameters
[in]xthe x values to use
[in]ythe function values to use
[in]nthe number of points to use to build the spline

The requested number of points must be at least the minimum required by the chosen interpolation method, and less than the minimum of the number of x values or number of y values. The first n points are then used to construct the interpolant. The resulting interpolant can subsequently be evaluated by calls to evaluate_at , evaluate_deriv_at , and evaluate_second_deriv_at .

Definition at line 255 of file gsl_interpolator.cpp.

References interp_type, interpolant, reset(), and to_gsl_interp_type().

double BubbleProfiler::GSL_interpolator::evaluate_at ( double  x) const

Evaluates computed interpolant at a point.

Parameters
[in]xthe point to evaluate at
Returns
the interpolated function value

The interpolant must previously have been computed via a call to construct_interpolant .

Definition at line 290 of file gsl_interpolator.cpp.

References accelerator, check_interpolant_computed(), check_within_bounds(), interpolant, and BubbleProfiler::make_raii_guard().

Referenced by set_interpolation_type().

double BubbleProfiler::GSL_interpolator::evaluate_deriv_at ( double  x) const

Evaluates derivative of computed interpolant at a point.

Parameters
[in]xthe point to evaluate at
Returns
the interpolated function derivative value

The interpolant must previously have been computed via a call to construct_interpolant .

Definition at line 317 of file gsl_interpolator.cpp.

References accelerator, check_interpolant_computed(), check_within_bounds(), interpolant, and BubbleProfiler::make_raii_guard().

Referenced by set_interpolation_type().

double BubbleProfiler::GSL_interpolator::evaluate_second_deriv_at ( double  x) const

Evaluates second derivative of computed interpolant at a point.

Parameters
[in]xthe point to evaluate at
Returns
the interpolated function second derivative value

The interpolant must previously have been computed via a call to construct_interpolant .

Definition at line 343 of file gsl_interpolator.cpp.

References accelerator, check_interpolant_computed(), check_within_bounds(), interpolant, and BubbleProfiler::make_raii_guard().

Referenced by set_interpolation_type().

GSL_interpolator & BubbleProfiler::GSL_interpolator::operator= ( const GSL_interpolator other)

Definition at line 202 of file gsl_interpolator.cpp.

References accelerator, interp_type, and interpolant.

GSL_interpolator& BubbleProfiler::GSL_interpolator::operator= ( GSL_interpolator &&  )
default
void BubbleProfiler::GSL_interpolator::reset ( Interpolation_type  t,
int  n_points 
)
private

Definition at line 402 of file gsl_interpolator.cpp.

References accelerator, interpolant, and to_gsl_interp_type().

Referenced by construct_interpolant().

void BubbleProfiler::GSL_interpolator::set_interpolation_type ( Interpolation_type  t)
inline

Sets the method used to perform interpolation.

Parameters
tthe interpolation method to use

Definition at line 91 of file gsl_interpolator.hpp.

References construct_interpolant(), evaluate_at(), evaluate_deriv_at(), evaluate_second_deriv_at(), and interp_type.

const gsl_interp_type * BubbleProfiler::GSL_interpolator::to_gsl_interp_type ( Interpolation_type  t) const
private

Member Data Documentation

std::unique_ptr<gsl_interp_accel, GSL_accel_deleter> BubbleProfiler::GSL_interpolator::accelerator {nullptr}
private
Interpolation_type BubbleProfiler::GSL_interpolator::interp_type {Interpolation_type::GSL_cspline}
private
std::unique_ptr<gsl_spline, GSL_spline_deleter> BubbleProfiler::GSL_interpolator::interpolant {nullptr}
private

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