|
NuMFor cf0f85d (2025-09-27)
Numerical (Modern) Fortran. Library for Simple Numerical computing
|
Type used to keep all information on splines. More...
Public Member Functions | |
| procedure, pass(csp) | cspl_interp (xc, csp) |
| Interpolates a function using previously calculated representation of splines. | |
| procedure, pass(csp) | cspl_interp_tab (xnew, csp) |
| Interpolates a function using previously calculated representation of splines. Works over an array of x values and returns an array of interpolated values | |
| generic | evaluate (xc, csp) |
| Interpolates a function using previously calculated representation of splines. | |
| generic | evaluate (xnew, csp) |
| Interpolates a function using previously calculated representation of splines. Works over an array of x values and returns an array of interpolated values | |
| procedure, pass(csp) | cspl_interpdev (xc, csp, m) |
| Interpolates the first derivative of a function. | |
| procedure, pass(csp) | cspl_interpdev_tab (xnew, csp, m) |
| generic | derivative (xc, csp, m) |
| Interpolates the first derivative of a function. | |
| generic | derivative (xnew, csp, m) |
| procedure, pass(csp) | integrate (xl, xu, csp, extrapolate) |
| Definite integral of a cubic spline function. | |
| procedure, pass(csp) | roots (csp) |
| csplroots Computes the roots of the Spline approximation | |
| type(cubicspline) function | init (x, y, s1, sn) |
Public Attributes | |
| real(dp), dimension(:, :), allocatable | s |
Coefficients of the polynomial S(:,i) valid for each interval i, x(i) <= x < x(i+1) | |
| real(dp), dimension(:), allocatable | x |
| Limits of the intervals. | |
Type used to keep all information on splines.
CubicSpline is the OO interface to cubic splines interpolation.
For each interval, the value will be 

And to evaluate the derivative
| procedure, pass(csp) cspl_interp | ( | real(dp), intent(in) | xc, |
| class(cubicspline), intent(in) | csp ) |
Interpolates a function using previously calculated representation of splines.
csplrep() | [in] | xc | value where evaluate the interpolation |
| [in] | csp | spline coefficients |
| procedure, pass(csp) cspl_interp_tab | ( | real(dp), dimension(:), intent(in) | xnew, |
| class(cubicspline), intent(in) | csp ) |
Interpolates a function using previously calculated representation of splines. Works over an array of x values and returns an array of interpolated values
| [in] | csp | spline coefficients |
| [in] | xnew | array of x values |
| procedure, pass(csp) cspl_interpdev | ( | real(dp), intent(in) | xc, |
| class(cubicspline), intent(in) | csp, | ||
| integer, intent(in), optional | m ) |
Interpolates the first derivative of a function.
csplrep() | [in] | xc | value where evaluate the interpolation |
| [in] | csp | spline coefficients |
| [in] | m | order of derivation |
| procedure, pass(csp) cspl_interpdev_tab | ( | real(dp), dimension(:), intent(in) | xnew, |
| class(cubicspline), intent(in) | csp, | ||
| integer, intent(in), optional | m ) |
| [in] | csp | spline coefficients |
| [in] | xnew | array of x values |
| [in] | m | order of derivative |
| generic derivative | ( | real(dp), intent(in) | xc, |
| class(cubicspline), intent(in) | csp, | ||
| integer, intent(in), optional | m ) |
Interpolates the first derivative of a function.
csplrep() | [in] | xc | value where evaluate the interpolation |
| [in] | csp | spline coefficients |
| [in] | m | order of derivation |
| generic derivative | ( | real(dp), dimension(:), intent(in) | xnew, |
| class(cubicspline), intent(in) | csp, | ||
| integer, intent(in), optional | m ) |
| [in] | csp | spline coefficients |
| [in] | xnew | array of x values |
| [in] | m | order of derivative |
| generic evaluate | ( | real(dp), intent(in) | xc, |
| class(cubicspline), intent(in) | csp ) |
Interpolates a function using previously calculated representation of splines.
csplrep() | [in] | xc | value where evaluate the interpolation |
| [in] | csp | spline coefficients |
| generic evaluate | ( | real(dp), dimension(:), intent(in) | xnew, |
| class(cubicspline), intent(in) | csp ) |
Interpolates a function using previously calculated representation of splines. Works over an array of x values and returns an array of interpolated values
| [in] | csp | spline coefficients |
| [in] | xnew | array of x values |
| type(cubicspline) function init | ( | real(dp), dimension(:), intent(in) | x, |
| real(dp), dimension(:), intent(in) | y, | ||
| real(dp), intent(in) | s1, | ||
| real(dp), intent(in) | sn ) |
| [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) |
References csplines::csplrep().

| procedure, pass(csp) integrate | ( | 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 |
| procedure, pass(csp) roots | ( | class(cubicspline), intent(in) | csp | ) |
csplroots Computes the roots of the Spline approximation
| [in] | csp | Spline approximation to consider |