NuMFor 9f2ab49 (2024-04-08)
Numerical (Modern) Fortran. Library for Simple Numerical computing
|
csplines implements interpolation using cubic splines Description: Submodule interpolate More...
Data Types | |
interface | csplev |
csplev Performs a spline interpolation in a point or in a table More... | |
interface | csplevder |
csplev Performs a spline interpolation in a point or in a table More... | |
interface | cubicspline |
Type used to keep all information on splines. More... | |
Functions/Subroutines | |
subroutine, public | csplrep (x, y, s1, sn, csp) |
cubic spline interpolation between tabulated data After calling this function the result may be used to evaluate the function as:ynew= csplev(xnew, csp) | |
type(cubicspline) function, public | csplder (csp, m) |
csplder Computes the derivative of the cubic spline | |
type(cubicspline) function, public | csplantider (csp, m) |
Computes the antiderivative of the CubicSpline approximation. | |
subroutine, public | cspl_clean (r) |
Clean-up a spline representation. | |
subroutine, public | cspleps (x, y, err) |
cspleps Estimates the error produced by using a spline approximation | |
real(dp) function, public | csplint (xl, xu, csp, extrapolate) |
Definite integral of a cubic spline function. | |
real(dp) function, public | csplint_square (xl, xu, csp) |
Integral of the square of a function expressed as a cubic spline. | |
real(dp) function, dimension(:), allocatable, public | csplroots (csp) |
csplroots Computes the roots of the Spline approximation | |
csplines implements interpolation using cubic splines Description: Submodule interpolate
subroutine, public cspl_clean | ( | type(cubicspline), intent(inout) | r | ) |
Clean-up a spline representation.
[in,out] | r | CubicSpline |
References basic::print_msg().
type(cubicspline) function, public csplantider | ( | type(cubicspline), intent(in) | csp, |
integer, intent(in) | m ) |
Computes the antiderivative of the CubicSpline approximation.
The result is a CubicSpline (not exactly, it is a polynomial of order m+4)
[in] | csp | CubicSpline object holding the spline |
[in] | m | order of integration |
References polynomial::polyint().
type(cubicspline) function, public csplder | ( | type(cubicspline), intent(in) | csp, |
integer, intent(in) | m ) |
csplder Computes the derivative of the cubic spline
[in] | csp | Interpolating object |
[in] | m | order of derivation (must be 1 or 2) |
References polynomial::polyder(), and basic::print_msg().
subroutine, public cspleps | ( | real(dp), dimension(:), intent(in) | x, |
real(dp), dimension(size(x)), intent(in) | y, | ||
real(dp), dimension(size(x)), intent(out) | err ) |
cspleps Estimates the error produced by using a spline approximation
This subroutine estimates the error introduced by natural cubic spline interpolation in a table x(i),y(i) (i=1,...,n)
. the interpolation error in the vicinity of x(k) is approximated by the difference between y(k) and the value obtained from the spline that interpolates the table with the k-th point removed. err is the largest relative error along the table.
[in] | x | grid points |
[in] | y | value of function at grid points |
[out] | err | Vector with error estimates |
real(dp) function, public csplint | ( | real(dp), intent(in) | xl, |
real(dp), intent(in) | xu, | ||
class(cubicspline), intent(in) | csp, | ||
logical, intent(in), optional | extrapolate ) |
Definite integral of a cubic spline function.
[in] | xl | Lower limit in the integral. |
[in] | xu | Upper limit in the integral. |
[in] | csp | Interpolating object |
[in] | extrapolate | Flag signaling if we extrapolate outside interval |
real(dp) function, public csplint_square | ( | real(dp), intent(in) | xl, |
real(dp), intent(in) | xu, | ||
type(cubicspline), intent(in) | csp ) |
Integral of the square of a function expressed as a cubic spline.
[in] | xl | Lower limit in the integral. |
[in] | xu | Upper limit in the integral. |
[in] | csp | Interpolating object |
subroutine, public csplrep | ( | real(dp), dimension(:), intent(in) | x, |
real(dp), dimension(:), intent(in) | y, | ||
real(dp), intent(in) | s1, | ||
real(dp), intent(in) | sn, | ||
type(cubicspline), intent(out) | csp ) |
cubic spline interpolation between tabulated data After calling this function the result may be used to evaluate the function as:
ynew= csplev(xnew, csp)
REF.: M.J. Maron, 'Numerical Analysis: A Practical Approach', Macmillan Publ. Co., New York 1982.
[in] | x | independent grid points |
[in] | y | corresponding function values |
[in] | s1 | second derivative at x(1) |
[in] | sn | second derivative at x(n) (Natural spline: s1=sn=0) |
[out] | csp | Coefficients stored in |
References basic::print_msg().
Referenced by cubicspline::init().
real(dp) function, dimension(:), allocatable, public csplroots | ( | class(cubicspline), intent(in) | csp | ) |
csplroots Computes the roots of the Spline approximation
[in] | csp | Spline approximation to consider |