OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
GSLSpline.h
Go to the documentation of this file.
1//
2// Copyright (c) 2023, Paul Scherrer Institute, Villigen PSI, Switzerland
3// All rights reserved
4//
5// This file is part of OPAL.
6//
7// OPAL is free software: you can redistribute it and/or modify
8// it under the terms of the GNU General Public License as published by
9// the Free Software Foundation, either version 3 of the License, or
10// (at your option) any later version.
11//
12// You should have received a copy of the GNU General License
13// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
14//
15
16#ifndef OPALX_GSL_SPLINE_H
17#define OPALX_GSL_SPLINE_H
18
19#include "CubicSpline.h"
20#include "LinearSpline.h"
21
22// Compatibility layer for GSL-like interface
23
25constexpr int gsl_interp_linear = 1;
26constexpr int gsl_interp_cspline = 0;
27
30
33
38inline gsl_spline* gsl_spline_alloc(const int type, size_t /*size*/) {
39 gsl_spline* result{};
40 if (type == gsl_interp_cspline) {
41 result = new CubicSpline();
42 } else if (type == gsl_interp_linear) {
43 result = new LinearSpline();
44 }
45 return result;
46}
47
51
57inline void gsl_spline_init(gsl_spline* spline, const double* x, const double* y, const size_t n) {
58 spline->init(std::vector<double>{x, x + n}, std::vector<double>{y, y + n});
59}
60
66inline double gsl_spline_eval(const gsl_spline* spline, const double x, gsl_interp_accel* accel) {
67 return spline->eval(x, *accel);
68}
69
77 const gsl_spline* spline, const double xa, const double xb, gsl_interp_accel* accel) {
78 return spline->evalIntegral(xa, xb, *accel);
79}
80
83inline void gsl_spline_free(const gsl_spline* spline) { delete spline; }
84
87inline void gsl_interp_accel_free(const gsl_interp_accel* accel) { delete accel; }
88
92 if (accel) {
93 accel->reset();
94 }
95}
96
97#endif // OPALX_GSL_SPLINE_H
double gsl_spline_eval_integ(const gsl_spline *spline, const double xa, const double xb, gsl_interp_accel *accel)
Evaluate the integral of a spline at x using an accelerator.
Definition GSLSpline.h:76
void gsl_interp_accel_reset(gsl_interp_accel *accel)
Reset an accelerator to the initial state.
Definition GSLSpline.h:91
constexpr int gsl_interp_cspline
Definition GSLSpline.h:26
double gsl_spline_eval(const gsl_spline *spline, const double x, gsl_interp_accel *accel)
Evaluate a spline at x using an accelerator.
Definition GSLSpline.h:66
gsl_interp_accel * gsl_interp_accel_alloc()
Allocate an interpolation accelerator.
Definition GSLSpline.h:50
void gsl_spline_init(gsl_spline *spline, const double *x, const double *y, const size_t n)
Initialize a spline with tabulated data.
Definition GSLSpline.h:57
void gsl_spline_free(const gsl_spline *spline)
Free a spline instance.
Definition GSLSpline.h:83
void gsl_interp_accel_free(const gsl_interp_accel *accel)
Free an accelerator instance.
Definition GSLSpline.h:87
gsl_spline * gsl_spline_alloc(const int type, size_t)
Allocate a spline instance (size ignored).
Definition GSLSpline.h:38
constexpr int gsl_interp_linear
GSL linear interpolation type identifiers.
Definition GSLSpline.h:25
Accelerator caching last interval indices.
Base class for the linear and cubic interpolation spline classes.
virtual double eval(double x, Accelerator &accel) const =0
Evaluate the spline at x using an accelerator.
virtual void init(const std::vector< double > &x, const std::vector< double > &y)
Initialize from tabulated data (natural spline).
virtual double evalIntegral(double xa, double xb, Accelerator &accel) const =0
Evaluate the integral of the spline at x.
Natural cubic spline interpolator.
Definition CubicSpline.h:28