fastmath.interpolation
1d, 2d interpolation functions.
See more:
Note: Smile interpolators doesn’t check ranges.
Input data
You provide data as sequence or double array.
1d interpolation
You provide two sequences:
xs
- x axis coorditanes, strictly monotonic (increasing)ys
- function values
See kriging-spline-interpolator
2d interpolation
This is grid based interpolation.
xs
- x axis coordinates, strictly monotonic (increasing)ys
- y axis coordinates, strictly monotonic (increasing)vs
- sequence of sequences of values (2d array) for all possible pairs. Array is column-wise:[ [first column] [second column] ...]
.
Categories
- Apache Commons Math interpolators: akima-spline bicubic divided-difference interpolators-1d-list interpolators-2d-list linear loess microsphere-2d-projection microsphere-projection neville piecewise-bicubic spline step-after step-before
- 1d interpolation: akima-spline cubic-spline divided-difference interpolators-1d-list kriging-spline linear linear-smile loess microsphere-projection neville rbf shepard spline step-after step-before
- 2d interpolation (grid based): bicubic bicubic-smile bilinear cubic-2d interpolators-2d-list microsphere-2d-projection piecewise-bicubic
- Smile interpolators: bicubic-smile bilinear cubic-2d cubic-spline interpolators-1d-list interpolators-2d-list kriging-spline linear-smile rbf shepard
Code snippets
Save graph
(defn save-graph
[f params & opts]
(let [fname (str "images/i/" (first opts) ".png")
[x1 x2] params]
(incanter.core/save (c/function-plot f x1 x2 :y-label "value")
(str "docs/" fname)
:width 250
:height 250)
fname))
Save interpolation graph
(defn save-interpolation
[f params & opts]
(let [xs [0.69 1.73 2.0 2.28 3.46 4.18 4.84 5.18 5.53 5.87 6.22]
ft (fn [x] (m/sin (* x (* 0.5 (m/cos (inc x))))))
ys (map ft xs)
[x1 x2] params
fname (str "images/i/" (first opts) ".png")]
(incanter.core/save (doto (c/function-plot ft 0 7 :y-label "value")
(c/set-theme-bw)
(c/add-points xs ys)
(c/add-function (f xs ys) x1 x2))
(str "docs/" fname)
:width 600
:height 300)
fname))
akima-spline
(akima-spline xs ys)
Create cubic spline interpolator using Akima algorithm. Minimum number of points: 5
xs[n] < xs[n+1] for all n.
Source: Apache Commons Math.
Examples
Akima plot
(save-interpolation akima-spline 0.69 6.22 ...)

cubic-2d
(cubic-2d xs ys vs)
Cubic spline 2d.
Grid based.
Source: Smile.
Examples
Usage
(let [interpolator
(cubic-2d [2 5 9] [2 3 10] [[4 0 2] [-1 2 -2] [-2 0 1]])]
(m/approx (interpolator 5.0 5.0)))
;;=> 4.68
;; Test: ok.
Array layout
(let [intrp (cubic-2d [2 5] [1 6] [[-1 -2] [3 4]])]
[(intrp 2 1) (intrp 2 6) (intrp 5 1) (intrp 5 6)])
;;=> [-1.0 -2.0 3.0 4.0]
cubic-spline
(cubic-spline xs ys)
Cubic spline interpolation.
Source: Smile.
Examples
Cubic spline plot
(save-interpolation cubic-spline 0 7 ...)

divided-difference
(divided-difference xs ys)
Create Divided Difference Algorithm for interpolation.
Source: Apache Commons Math.
Examples
Divided Difference plot
(save-interpolation divided-difference 0.65 6.5 ...)

interpolators-1d-list
Map of 1d interpolation functions
Examples
List of names
(keys interpolators-1d-list)
;;=> (:linear-smile :divided-difference
;;=> :rbf :neville
;;=> :step-after :kriging-spline
;;=> :cubic-spline :spline
;;=> :loess :step-before
;;=> :shepard :linear
;;=> :microsphere :akima)
interpolators-2d-list
Map of 2d interpolation functions
Examples
List of names
(keys interpolators-2d-list)
;;=> (:bicubic :piecewise-bicubic
;;=> :microsphere-2d :bilinear
;;=> :bicubic-smile :cubic-2d)
kriging-spline
(kriging-spline xs ys)
Kriging interpolation.
Source: Smile.
Examples
Usage
(let [interpolator (kriging-spline [2 5 9 10 11]
[0.4 1.0 -1.0 -0.5 0.0])]
(m/approx (interpolator 7.0)))
;;=> -0.07
;; Test: ok.
Kriging spline plot
(save-interpolation kriging-spline 0 7 ...)

linear
(linear xs ys)
Create Divided Difference Algorithm for inqterpolation.
Source: Apache Commons Math.
Examples
Linear (Apache Commons version) plot
(save-interpolation linear 0.69 6.22 ...)

linear-smile
(linear-smile xs ys)
Linear interpolation from Smile library.
Source: Smile.
Examples
Linear (SMILE version) plot
(save-interpolation linear-smile 0 7 ...)

loess
(loess xs ys)
(loess bandwidth robustness-iters xs ys)
(loess bandwidth robustness-iters accuracy xs ys)
Local Regression Algorithm
- bandwidth: 0.2-1.0 (optimal: 0.25-0.5, default: 0.4)
- robustness-iters: 0-4 (optimal: 0, default: 2)
- accuracy: double (default: 1e-12)
Source: Apache Commons Math.
Examples
Loess plot
(save-interpolation loess 0.69 6.22 ...)

Loess plot, bandwidth=0.7, robustness-iters=4
(save-interpolation (fn [xs ys] (loess 0.7 4 xs ys)) 0.69 6.22 ...)

Loess plot, bandwidth=0.2, robustness-iters=1
(save-interpolation (fn [xs ys] (loess 0.2 1 xs ys)) 0.69 6.22 ...)

microsphere-2d-projection
(microsphere-2d-projection elements max-dark-friction dark-threshold background exponent shared-sphere? no-interpolation-tolerance xs ys vs)
Microsphere projection interpolator - 2d version
Grid based.
Source: Apache Commons Math.
microsphere-projection
(microsphere-projection elements max-dark-friction dark-threshold background exponent shared-sphere? no-interpolation-tolerance xs ys)
Microsphere projection interpolator - 1d version
Source: Apache Commons Math.
Examples
Microsphere projection plot
(save-interpolation
(fn [xs ys]
(microsphere-projection 4 0.5 1.0E-7 0.1 2.5 false 0.1 xs ys))
0
7
...)

neville
(neville xs ys)
Neville algorithm
Source: Apache Commons Math.
Examples
Neville plot
(save-interpolation neville 0.65 6.5 ...)

piecewise-bicubic
(piecewise-bicubic xs ys vs)
Piecewise bicubic 2d.
Grid based.
Source: Apache Commons Math.
rbf
(rbf rbf-fn normalize? xs ys)
(rbf rbf-fn xs ys)
RBF (Radial Basis Function) interpolation.
Source: Smile
Examples
RBF - Linear
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :linear) xs ys)) 0 7 ...)

RBF - Gaussian
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :gaussian) xs ys))
0
7
...)

RBF - Multiquadratic
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :multiquadratic) xs ys))
0
7
...)

RBF - Inverse ultiquadratic
(save-interpolation (fn [xs ys]
(rbf (rbf/rbf :inverse-multiquadratic) xs ys))
0
7
...)

RBF - Inverse quadratic
(save-interpolation (fn [xs ys]
(rbf (rbf/rbf :inverse-quadratic) xs ys))
0
7
...)

RBF - Thinplate
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :thinplate) xs ys))
0
7
...)

RBF - Polyharmonic 3
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :polyharmonic 3) xs ys))
0
7
...)

RBF - Polyharmonic 4
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :polyharmonic 4) xs ys))
0
7
...)

RBF - Wendland
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :wendland) xs ys))
0
7
...)

RBF - Wu
(save-interpolation (fn [xs ys] (rbf (rbf/rbf :wu) xs ys)) 0 7 ...)

shepard
(shepard xs ys)
(shepard p xs ys)
Shepard interpolation.
Source: Smile.
Examples
Shepard plot
(save-interpolation shepard 0 7 ...)

Shepard plot, p=5
(save-interpolation (fn [xs ys] (shepard 5 xs ys)) 0 7 ...)

Shepard plot, p=0.9
(save-interpolation (fn [xs ys] (shepard 0.9 xs ys)) 0 7 ...)

spline
(spline xs ys)
Cubic spline interpolation
Source: Apache Commons Math.
Examples
Spline plot
(save-interpolation spline 0.69 6.22 ...)

step-after
(step-after xs ys)
Step function.
Examples
Step function plot
(save-interpolation step-after 0 7 ...)

step-before
(step-before xs ys)
Step function.
Examples
Step function plot
(save-interpolation step-before 0 7 ...)
