Reference card for the numeric module
|
Function | Description
|
| abs | Absolute value
| acos | Arc-cosine
| add | Pointwise sum x+y
| addeq | Pointwise sum x+=y
| all | All the components of x are true
| and | Pointwise x && y
| andeq | Pointwise x &= y
| any | One or more of the components of x are true
| asin | Arc-sine
| atan | Arc-tangeant
| atan2 | Arc-tangeant (two parameters)
| band | Pointwise x & y
| bench | Benchmarking routine
| bnot | Binary negation ~x
| bor | Binary or x|y
| bxor | Binary xor x^y
| cLU | Coordinate matrix LU decomposition
| cLUsolve | Coordinate matrix LU solve
| cdelsq | Coordinate matrix Laplacian
| cdotMV | Coordinate matrix-vector product
| ceil | Pointwise Math.ceil(x)
| cgrid | Coordinate grid for cdelsq
| clone | Deep copy of Array
| cos | Pointwise Math.cos(x)
| det | Determinant
| diag | Create diagonal matrix
| dim | Get Array dimensions
| div | Pointwise x/y
| diveq | Pointwise x/=y
| dopri | Numerical integration of ODE using Dormand-Prince RK method. Returns an object Dopri.
| Dopri.at | Evaluate the ODE solution at a point
| dot | Matrix-Matrix, Matrix-Vector and Vector-Matrix product
| eig | Eigenvalues and eigenvectors
| epsilon | 2.220446049250313e-16
| eq | Pointwise comparison x === y
| exp | Pointwise Math.exp(x)
| floor | Poinwise Math.floor(x)
| geq | Pointwise x>=y
| getBlock | Extract a block from a matrix
| getDiag | Get the diagonal of a matrix
| gt | Pointwise x>y
| identity | Identity matrix
| imageURL | Encode a matrix as an image URL
|
|
Function | Description
|
| inv | Matrix inverse
| isFinite | Pointwise isFinite(x)
| isNaN | Pointwise isNaN(x)
| largeArray | Don't prettyPrint Arrays larger than this
| leq | Pointwise x<=y
| linspace | Generate evenly spaced values
| log | Pointwise Math.log(x)
| lshift | Pointwise x<<y
| lshifteq | Pointwise x<<=y
| lt | Pointwise x<y
| mapreduce | Make a pointwise map-reduce function
| mod | Pointwise x%y
| modeq | Pointwise x%=y
| mul | Pointwise x*y
| neg | Pointwise -x
| neq | Pointwise x!==y
| norm2 | Square root of the sum of the square of the entries of x
| norm2Squared | Sum of squares of entries of x
| norminf | Largest modulus entry of x
| not | Pointwise logical negation !x
| or | Pointwise logical or x||y
| oreq | Pointwise x|=y
| parseCSV | Parse a CSV file into an Array
| parseDate | Pointwise parseDate(x)
| parseFloat | Pointwise parseFloat(x)
| pointwise | Create a pointwise function
| pow | Pointwise Math.pow(x)
| precision | Number of digits to prettyPrint
| prettyPrint | Pretty-prints x
| random | Create an Array of random numbers
| rep | Create an Array by duplicating values
| round | Pointwise Math.round(x)
| rrshift | Pointwise x>>>y
| rrshifteq | Pointwise x>>>=y
| rshift | Pointwise x>>y
| rshifteq | Pointwise x>>=y
| sLUP | Sparse LUP decomposition
| sLUPsolve | Sparse LUP solve
| same | x and y are entrywise identical
| sclone | Sparse clone Array
| sdiag | Sparse diagonal matrix
| sdim | Sparse dimension of Array
|
|
Function | Description
|
| sdot | Sparse matrix-matrix, matrix-vector, vector-matrix and vector-vector product.
| seedrandom | The seedrandom module
| setBlock | Set a block of a matrix
| sgather | Sparse gather values
| sidentity | Sparse identity
| sin | Pointwise Math.sin(x)
| solveQP | Solve a quadratic programming problem
| spline | Create a Spline object
| Spline.at | Evaluate the Spline at a point
| Spline.diff | Differentiate the Spline
| Spline.roots | Find all the roots of the Spline
| sqrt | Pointwise Math.sqrt(x)
| sscatter | Sparse scatter values
| stranspose | Sparse transpose
| sub | Pointwise x-y
| subeq | Pointwise x-=y
| sum | Sum all the entries of x
| svd | Singular value decomposition
| t | Create a tensor type T (may be complex-valued)
| T.<numericfun> | Supported <numericfun> are: abs, add, cos, diag, div, dot, exp, getBlock, getDiag, inv, log, mul, neg, norm2, setBlock, sin, sub, transpose
| T.conj | Pointwise complex conjugate
| T.fft | Fast Fourier transform
| T.get | Read an entry
| T.getRow | Get a row
| T.getRows | Get a range of rows
| T.ifft | Inverse FFT
| T.reciprocal | Pointwise 1/z
| T.set | Set an entry
| T.setRow | Set a row
| T.setRows | Set a range of rows
| T.transjugate | The conjugate-transpose of a matrix
| tan | Pointwise Math.tan(x)
| tensor | Tensor product ret[i][j] = x[i]*y[j]
| toCSV | Make a CSV file
| transpose | Matrix transpose
| uncmin | Unconstrained optimization
| version | Version string for the numeric library
| xor | Pointwise x^y
| xoreq | Pointwise x^=y
| |
Numerical analysis in Javascript
Numeric Javascript is
library that provides many useful functions for numerical
calculations, particularly for linear algebra (vectors and matrices).
You can create vectors and matrices and multiply them:
> A = [[1,2,3],[4,5,6]];
[[1,2,3],
[4,5,6]]
> x = [7,8,9]
[7,8,9]
> numeric.dot(A,x);
[50,122]
The example shown above can be executed in the
Javascript Workshop or at any
Javascript prompt. The Workshop provides plotting capabilities:

The function workshop.plot() is essentially the flot
plotting command.
The numeric library provides functions that implement most of the usual Javascript
operators for vectors and matrices:
> y = [10,1,2];
[10,1,2]
> numeric['+'](x,y)
[17,9,11]
> numeric['>'](x,y)
[false,true,true]
These operators can also be called with plain Javascript function names:
> numeric.add(x,y)
[17,9,11]
You can also use these operators with three or more parameters:
> numeric.add([1,2],[3,4],[5,6],[7,8])
[16,20]
The function numeric.inv() can be used to compute the inverse of an invertible matrix:
> A = [[1,2,3],[4,5,6],[7,1,9]]
[[1,2,3],
[4,5,6],
[7,1,9]]
> Ainv = numeric.inv(A);
[[-0.9286,0.3571,0.07143],
[-0.1429,0.2857,-0.1429],
[0.7381,-0.3095,0.07143]]
The function numeric.prettyPrint() is used to print most of the examples in this documentation.
It formats objects, arrays and numbers so that they can be understood easily. All output is automatically
formatted using numeric.prettyPrint() when in the
Workshop. In order to present the information clearly and
succintly, the function numeric.prettyPrint() lays out matrices so that all the numbers align.
Furthermore, numbers are given approximately using the numeric.precision variable:
> numeric.precision = 10; x = 3.141592653589793
3.141592654
> numeric.precision = 4; x
3.142
The default precision is 4 digits. In addition to printing approximate numbers,
the function numeric.prettyPrint() will replace large arrays with the string ...Large Array...:
> numeric.identity(100)
...Large Array...
By default, this happens with the Array's length is more than 50. This can be controlled by setting the
variable numeric.largeArray to an appropriate value:
> numeric.largeArray = 2; A = numeric.identity(4)
...Large Array...
> numeric.largeArray = 50; A
[[1,0,0,0],
[0,1,0,0],
[0,0,1,0],
[0,0,0,1]]
In particular, if you want to print all Arrays regardless of size, set numeric.largeArray = Infinity.
Math Object functions
The Math object functions have also been adapted to work on Arrays as follows:
> numeric.exp([1,2]);
[2.718,7.389]
> numeric.exp([[1,2],[3,4]])
[[2.718, 7.389],
[20.09, 54.6]]
> numeric.abs([-2,3])
[2,3]
> numeric.acos([0.1,0.2])
[1.471,1.369]
> numeric.asin([0.1,0.2])
[0.1002,0.2014]
> numeric.atan([1,2])
[0.7854,1.107]
> numeric.atan2([1,2],[3,4])
[0.3218,0.4636]
> numeric.ceil([-2.2,3.3])
[-2,4]
> numeric.floor([-2.2,3.3])
[-3,3]
> numeric.log([1,2])
[0,0.6931]
> numeric.pow([2,3],[0.25,7.1])
[1.189,2441]
> numeric.round([-2.2,3.3])
[-2,3]
> numeric.sin([1,2])
[0.8415,0.9093]
> numeric.sqrt([1,2])
[1,1.414]
> numeric.tan([1,2])
[1.557,-2.185]
Utility functions
The function numeric.dim() allows you to compute the dimensions of an Array.
> numeric.dim([1,2])
[2]
> numeric.dim([[1,2,3],[4,5,6]])
[2,3]
You can perform a deep comparison of Arrays using numeric.same():
> numeric.same([1,2],[1,2])
true
> numeric.same([1,2],[1,2,3])
false
> numeric.same([1,2],[[1],[2]])
false
> numeric.same([[1,2],[3,4]],[[1,2],[3,4]])
true
> numeric.same([[1,2],[3,4]],[[1,2],[3,5]])
false
> numeric.same([[1,2],[2,4]],[[1,2],[3,4]])
false
You can create a multidimensional Array from a given value using numeric.rep()
> numeric.rep([3],5)
[5,5,5]
> numeric.rep([2,3],0)
[[0,0,0],
[0,0,0]]
You can loop over Arrays as you normally would. However, in order to quickly generate optimized
loops, the numeric library provides a few efficient loop-generation mechanisms. For example, the
numeric.mapreduce() function can be used to make a function that computes the sum of all the
entries of an Array.
> sum = numeric.mapreduce('accum += xi','0'); sum([1,2,3])
6
> sum([[1,2,3],[4,5,6]])
21
The functions numeric.any() and numeric.all() allow you to check whether any or all entries
of an Array are boolean true values.
> numeric.any([false,true])
true
> numeric.any([[0,0,3.14],[0,false,0]])
true
> numeric.any([0,0,false])
false
> numeric.all([false,true])
false
> numeric.all([[1,4,3.14],["no",true,-1]])
true
> numeric.all([0,0,false])
false
You can also create "pointwise" functions. This is how you create a vector sum function:
> add = numeric.pointwise(['x[i]','y[i]'],'ret[i] = x[i]+y[i];'); add([1,2],[3,4])
[4,6]
You can create a diagonal matrix using numeric.diag()
> numeric.diag([1,2,3])
[[1,0,0],
[0,2,0],
[0,0,3]]
The function numeric.identity() returns the identity matrix.
> numeric.identity(3)
[[1,0,0],
[0,1,0],
[0,0,1]]
Random Arrays can also be created:
> numeric.random([2,3])
[[0.05303,0.1537,0.7280],
[0.3839,0.08818,0.6316]]
You can generate a vector of evenly spaced values:
> numeric.linspace(1,5);
[1,2,3,4,5]
> numeric.linspace(1,3,5);
[1,1.5,2,2.5,3]
Arithmetic operations
The standard arithmetic operations have been vectorized:
> numeric.addVV([1,2],[3,4])
[4,6]
> numeric.addVS([1,2],3)
[4,5]
There are also polymorphic functions:
> numeric.add(1,[2,3])
[3,4]
> numeric.add([1,2,3],[4,5,6])
[5,7,9]
The other arithmetic operations are available:
> numeric.sub([1,2],[3,4])
[-2,-2]
> numeric.mul([1,2],[3,4])
[3,8]
> numeric.div([1,2],[3,4])
[0.3333,0.5]
The in-place operators (such as +=) are also available:
> v = [1,2,3,4]; numeric.addeq(v,3); v
[4,5,6,7]
> numeric.subeq([1,2,3],[5,3,1])
[-4,-1,2]
Unary operators:
> numeric.neg([1,-2,3])
[-1,2,-3]
> numeric.isFinite([10,NaN,Infinity])
[true,false,false]
> numeric.isNaN([10,NaN,Infinity])
[false,true,false]
Linear algebra
Matrix products are implemented in the functions
numeric.dotVV()
numeric.dotVM()
numeric.dotMV()
numeric.dotMM():
> numeric.dotVV([1,2],[3,4])
11
> numeric.dotVM([1,2],[[3,4],[5,6]])
[13,16]
> numeric.dotMV([[1,2],[3,4]],[5,6])
[17,39]
> numeric.dotMMbig([[1,2],[3,4]],[[5,6],[7,8]])
[[19,22],
[43,50]]
> numeric.dotMMsmall([[1,2],[3,4]],[[5,6],[7,8]])
[[19,22],
[43,50]]
> numeric.dot([1,2,3,4,5,6,7,8,9],[1,2,3,4,5,6,7,8,9])
285
The function numeric.dot() is "polymorphic" and selects the appropriate Matrix product:
> numeric.dot([1,2,3],[4,5,6])
32
> numeric.dot([[1,2,3],[4,5,6]],[7,8,9])
[50,122]
The determinant:
> numeric.det([[1,2],[3,4]]);
-2
> numeric.det([[6,8,4,2,8,5],[3,5,2,4,9,2],[7,6,8,3,4,5],[5,5,2,8,1,6],[3,2,2,4,2,2],[8,3,2,2,4,1]]);
-1404
The matrix inverse:
> numeric.inv([[1,2],[3,4]])
[[ -2, 1],
[ 1.5, -0.5]]
The transpose:
> numeric.transpose([[1,2,3],[4,5,6]])
[[1,4],
[2,5],
[3,6]]
> numeric.transpose([[1,2,3,4,5,6,7,8,9,10,11,12]])
[[ 1],
[ 2],
[ 3],
[ 4],
[ 5],
[ 6],
[ 7],
[ 8],
[ 9],
[10],
[11],
[12]]
You can compute the 2-norm of an Array, which is the square root of the sum of the squares of the entries.
> numeric.norm2([1,2])
2.236
Computing the tensor product of two vectors:
> numeric.tensor([1,2],[3,4,5])
[[3,4,5],
[6,8,10]]
Data manipulation
There are also some data manipulation functions. You can parse dates:
> numeric.parseDate(['1/13/2013','2001-5-9, 9:31']);
[1.358e12,9.894e11]
Parse floating point quantities:
> numeric.parseFloat(['12','0.1'])
[12,0.1]
Parse CSV files:
> numeric.parseCSV('a,b,c\n1,2.3,.3\n4e6,-5.3e-8,6.28e+4')
[[ "a", "b", "c"],
[ 1, 2.3, 0.3],
[ 4e6, -5.3e-8, 62800]]
> numeric.toCSV([[1.23456789123,2],[3,4]])
"1.23456789123,2
3,4
"
You can also fetch a URL (a thin wrapper around XMLHttpRequest):
> numeric.getURL('tools/helloworld.txt').responseText
"Hello, world!"
Complex linear algebra
You can also manipulate complex numbers:
> z = new numeric.T(3,4);
{x: 3, y: 4}
> z.add(5)
{x: 8, y: 4}
> w = new numeric.T(2,8);
{x: 2, y: 8}
> z.add(w)
{x: 5, y: 12}
> z.mul(w)
{x: -26, y: 32}
> z.div(w)
{x:0.5588,y:-0.2353}
> z.sub(w)
{x:1, y:-4}
Complex vectors and matrices can also be handled:
> z = new numeric.T([1,2],[3,4]);
{x: [1,2], y: [3,4]}
> z.abs()
{x:[3.162,4.472],y:}
> z.conj()
{x:[1,2],y:[-3,-4]}
> z.norm2()
5.477
> z.exp()
{x:[-2.691,-4.83],y:[0.3836,-5.592]}
> z.cos()
{x:[-1.528,-2.459],y:[0.1658,-2.745]}
> z.sin()
{x:[0.2178,-2.847],y:[1.163,2.371]}
> z.log()
{x:[1.151,1.498],y:[1.249,1.107]}
Complex matrices:
> A = new numeric.T([[1,2],[3,4]],[[0,1],[2,-1]]);
{x:[[1, 2],
[3, 4]],
y:[[0, 1],
[2,-1]]}
> A.inv();
{x:[[0.125,0.125],
[ 0.25, 0]],
y:[[ 0.5,-0.25],
[-0.375,0.125]]}
> A.inv().dot(A)
{x:[[1, 0],
[0, 1]],
y:[[0,-2.776e-17],
[0, 0]]}
> A.get([1,1])
{x: 4, y: -1}
> A.transpose()
{ x: [[1, 3],
[2, 4]],
y: [[0, 2],
[1,-1]] }
> A.transjugate()
{ x: [[ 1, 3],
[ 2, 4]],
y: [[ 0,-2],
[-1, 1]] }
> numeric.T.rep([2,2],new numeric.T(2,3));
{ x: [[2,2],
[2,2]],
y: [[3,3],
[3,3]] }
Eigenvalues
Eigenvalues:
> A = [[1,2,5],[3,5,0],[7,-3,5]];
[[ 1, 2, 5],
[ 3, 5, 0],
[ 7,-3, 5]]
> B = numeric.eig(A);
{lambda:{x:[-4.068,8.742,6.326],y:},
E: {x: [[ 0.7281,-0.5634, 0.4034 ],
[-0.2409,-0.4516, 0.9127 ],
[-0.6417,-0.6918, 0.06462]],
y:}}
> C = B.E.dot(numeric.T.diag(B.lambda)).dot(B.E.inv());
{x:[[ 1, 2, 5 ],
[ 3, 5,-3.109e-15],
[ 7, -3, 5 ]],
y:}
> C.sub(A).norm2()<1e-13
true
Note that eigenvalues and eigenvectors are returned as complex numbers (type numeric.T). This is because
eigenvalues are often complex even when the matrix is real.
Singular value decomposition (Shanti Rao)
Shanti Rao kindly emailed me an implementation of the Golub and Reisch algorithm:
> A = [[22,10,2,3,7],[14,7,10,0,8],[-1,13,-1,-11,3],[-3,-2,13,-2,4],[9,8,1,-2,4],[9,1,-7,5,-1],[2,-6,6,5,1],[4,5,0,-2,2]]
[[ 22, 10, 2, 3, 7],
[ 14, 7, 10, 0, 8],
[ -1, 13, -1,-11, 3],
[ -3, -2, 13, -2, 4],
[ 9, 8, 1, -2, 4],
[ 9, 1, -7, 5, -1],
[ 2, -6, 6, 5, 1],
[ 4, 5, 0, -2, 2]]
> numeric.svd(A)
{U:
[[ -0.7071, -0.1581, 0.1768, 0.2494, 0.4625],
[ -0.5303, -0.1581, -0.3536, 0.1556, -0.4984],
[ -0.1768, 0.7906, -0.1768, -0.1546, 0.3967],
[ -1.506e-17, -0.1581, -0.7071, -0.3277, 0.1],
[ -0.3536, 0.1581, 1.954e-15, -0.07265, -0.2084],
[ -0.1768, -0.1581, 0.5303, -0.5726, -0.05555],
[ -7.109e-18, -0.4743, -0.1768, -0.3142, 0.4959],
[ -0.1768, 0.1581, 1.915e-15, -0.592, -0.2791]],
S:
[ 35.33, 20, 19.6, 0, 0],
V:
[[ -0.8006, -0.3162, 0.2887, -0.4191, 0],
[ -0.4804, 0.6325, 7.768e-15, 0.4405, 0.4185],
[ -0.1601, -0.3162, -0.866, -0.052, 0.3488],
[ 4.684e-17, -0.6325, 0.2887, 0.6761, 0.2442],
[ -0.3203, 3.594e-15, -0.2887, 0.413, -0.8022]]}
Sparse linear algebra
Sparse linear algebra:
> numeric.sidentity(3)
[[1],
[ ,1],
[ , ,1]]
> numeric.stranspose([[1],[2,3],[4,5,6]])
[[1,2,4],
[ ,3,5],
[ , ,6]]
> A = [[2,-1],[-1,2,-1],[,-1,2]]; lup = numeric.sLUP(A)
{L: [[ 1],
[-0.5, 1],
[ ,-0.6667, 1]],
U: [[ 2, -1],
[ , 1.5, -1],
[ , , 1.333]],
P: [0,1,2],
Pinv: [0,1,2]}
> numeric.sdot(lup.L,lup.U)
[[ 2, -1],
[ -1, 2, -1],
[ , -1, 2]]
> x = [3,1,7]; b = numeric.sdot(A,x);
[5,-8,13]
> numeric.sLUPsolve(lup,b)
[3,1,7]
The numeric.sscatter() and numeric.sgather() functions can be used to convert between
sparse matrices and the coordinate encoding:
> A = numeric.sscatter([[0,0,1,1,1,2,2],[0,1,0,1,2,1,2],[1,2,3,4,5,6,7]])
[[1,2],
[3,4,5],
[ ,6,7]]
> numeric.sgather(A)
[[0,0,1,1,1,2,2],
[0,1,0,1,2,1,2],
[1,2,3,4,5,6,7]]
Coordinate matrices
Sparse matrices can be very slow in certain browsers. Functions that
work with matrices in the coordinate encoding using band strategies may execute faster if
your matrices happen to be banded.
LU decomposition:
> lu = numeric.cLU([[0,0,1,1,1,2,2],[0,1,0,1,2,1,2],[2,-1,-1,2,-1,-1,2]])
{U:[[ 0, 0, 1, 1, 2 ],
[ 0, 1, 1, 2, 2 ],
[ 2, -1, 1.5, -1, 1.333]],
L:[[ 0, 1, 1, 2, 2 ],
[ 0, 0, 1, 1, 2 ],
[ 1, -0.5, 1,-0.6667, 1 ]]}
> numeric.cLUsolve(lu,[5,-8,13])
[3,1,7]
Note that numeric.cLU() does not have any pivoting.
Solving PDEs
The functions numeric.cgrid() and numeric.cdelsq() can be used to obtain a
numerical Laplacian for a domain.
> g = numeric.cgrid(5)
[[-1,-1,-1,-1,-1],
[-1, 0, 1, 2,-1],
[-1, 3, 4, 5,-1],
[-1, 6, 7, 8,-1],
[-1,-1,-1,-1,-1]]
> coordL = numeric.cdelsq(g)
[[ 0, 0, 0, 1, 1, 1, 1, 2, 2, 2, 3, 3, 3, 3, 4, 4, 4, 4, 4, 5, 5, 5, 5, 6, 6, 6, 7, 7, 7, 7, 8, 8, 8],
[ 1, 3, 0, 0, 2, 4, 1, 1, 5, 2, 0, 4, 6, 3, 1, 3, 5, 7, 4, 2, 4, 8, 5, 3, 7, 6, 4, 6, 8, 7, 5, 7, 8],
[-1,-1, 4,-1,-1,-1, 4,-1,-1, 4,-1,-1,-1, 4,-1,-1,-1,-1, 4,-1,-1,-1, 4,-1,-1, 4,-1,-1,-1, 4,-1,-1, 4]]
> L = numeric.sscatter(coordL); // Just to see what it looks like
[[ 4, -1, , -1],
[ -1, 4, -1, , -1],
[ , -1, 4, , , -1],
[ -1, , , 4, -1, , -1],
[ , -1, , -1, 4, -1, , -1],
[ , , -1, , -1, 4, , , -1],
[ , , , -1, , , 4, -1],
[ , , , , -1, , -1, 4, -1],
[ , , , , , -1, , -1, 4]]
> lu = numeric.cLU(coordL); x = numeric.cLUsolve(lu,[1,1,1,1,1,1,1,1,1]);
[0.6875,0.875,0.6875,0.875,1.125,0.875,0.6875,0.875,0.6875]
> numeric.cdotMV(coordL,x)
[1,1,1,1,1,1,1,1,1]
> G = numeric.rep([5,5],0); for(i=0;i<5;i++) for(j=0;j<5;j++) if(g[i][j]>=0) G[i][j] = x[g[i][j]]; G
[[ 0 , 0 , 0 , 0 , 0 ],
[ 0 , 0.6875, 0.875 , 0.6875, 0 ],
[ 0 , 0.875 , 1.125 , 0.875 , 0 ],
[ 0 , 0.6875, 0.875 , 0.6875, 0 ],
[ 0 , 0 , 0 , 0 , 0 ]]
> workshop.html('<img src="'+numeric.imageURL(numeric.mul([G,G,G],200))+'" width=100 />');
You can also work on an L-shaped or arbitrary-shape domain:
> numeric.cgrid(6,'L')
[[-1,-1,-1,-1,-1,-1],
[-1, 0, 1,-1,-1,-1],
[-1, 2, 3,-1,-1,-1],
[-1, 4, 5, 6, 7,-1],
[-1, 8, 9,10,11,-1],
[-1,-1,-1,-1,-1,-1]]
> numeric.cgrid(5,function(i,j) { return i!==2 || j!==2; })
[[-1,-1,-1,-1,-1],
[-1, 0, 1, 2,-1],
[-1, 3,-1, 4,-1],
[-1, 5, 6, 7,-1],
[-1,-1,-1,-1,-1]]
Cubic splines
You can do some (natural) cubic spline interpolation:
> numeric.spline([1,2,3,4,5],[1,2,1,3,2]).at(numeric.linspace(1,5,10))
[ 1, 1.731, 2.039, 1.604, 1.019, 1.294, 2.364, 3.085, 2.82, 2]
For clamped splines:
> numeric.spline([1,2,3,4,5],[1,2,1,3,2],0,0).at(numeric.linspace(1,5,10))
[ 1, 1.435, 1.98, 1.669, 1.034, 1.316, 2.442, 3.017, 2.482, 2]
For periodic splines:
> numeric.spline([1,2,3,4],[0.8415,0.04718,-0.8887,0.8415],'periodic').at(numeric.linspace(1,4,10))
[ 0.8415, 0.9024, 0.5788, 0.04718, -0.5106, -0.8919, -0.8887, -0.3918, 0.3131, 0.8415]
Vector splines:
> numeric.spline([1,2,3],[[0,1],[1,0],[0,1]]).at(1.78)
[0.9327,0.06728]
Spline differentiation:
> xs = [0,1,2,3]; numeric.spline(xs,numeric.sin(xs)).diff().at(1.5)
0.07302
Find all the roots:
> xs = numeric.linspace(0,30,31); ys = numeric.sin(xs); s = numeric.spline(xs,ys).roots();
[0, 3.142, 6.284, 9.425, 12.57, 15.71, 18.85, 21.99, 25.13, 28.27]
Fast Fourier Transforms
FFT on numeric.T objects:
> z = (new numeric.T([1,2,3,4,5],[6,7,8,9,10])).fft()
{x:[15,-5.941,-3.312,-1.688, 0.941],
y:[40, 0.941,-1.688,-3.312,-5.941]}
Also inverse FFT:
> z.ifft()
{x:[1,2,3,4,5],
y:[6,7,8,9,10]}
Quadratic Programming (Alberto Santini)
The Quadratic Programming function numeric.solveQP() is based on Alberto Santini's
quadprog, which is itself a port of the corresponding
R routines.
> numeric.solveQP([[1,0,0],[0,1,0],[0,0,1]],[0,5,0],[[-4,2,0],[-3,1,-2],[0,0,1]],[-8,2,0]);
{ solution: [0.4762,1.048,2.095],
value: [-2.381],
unconstrained_solution:[ 0, 5, 0],
iterations: [ 3, 0],
iact: [ 3, 2, 0],
message: "" }
Unconstrained optimization
We also include a simple unconstrained optimization routine. Here are some demos from from Moré et al., 1981:
> sqr = function(x) { return x*x; }; numeric.uncmin(function(x) { return sqr(10*(x[1]-x[0]*x[0])) + sqr(1-x[0]); },[-1.2,1]).solution
[1,1]
> f = function(x) { return sqr(-13+x[0]+((5-x[1])*x[1]-2)*x[1])+sqr(-29+x[0]+((x[1]+1)*x[1]-14)*x[1]); }; x0 = numeric.uncmin(f,[0.5,-2]).solution
[11.41,-0.8968]
> f = function(x) { return sqr(1e4*x[0]*x[1]-1)+sqr(Math.exp(-x[0])+Math.exp(-x[1])-1.0001); }; x0 = numeric.uncmin(f,[0,1]).solution
[1.098e-5,9.106]
> f = function(x) { return sqr(x[0]-1e6)+sqr(x[1]-2e-6)+sqr(x[0]*x[1]-2)}; x0 = numeric.uncmin(f,[0,1]).solution
[1e6,2e-6]
> f = function(x) { return sqr(1.5-x[0]*(1-x[1]))+sqr(2.25-x[0]*(1-x[1]*x[1]))+sqr(2.625-x[0]*(1-x[1]*x[1]*x[1])); }; x0 = numeric.uncmin(f,[1,1]).solution
[3,0.5]
> f = function(x) { var ret = 0,i; for(i=1;i<=10;i++) ret+=sqr(2+2*i-Math.exp(i*x[0])-Math.exp(i*x[1])); return ret; }; x0 = numeric.uncmin(f,[0.3,0.4]).solution
[0.2578,0.2578]
> y = [0.14,0.18,0.22,0.25,0.29,0.32,0.35,0.39,0.37,0.58,0.73,0.96,1.34,2.10,4.39]; f = function(x) { var ret = 0,i; for(i=1;i<=15;i++) ret+=sqr(y[i-1]-(x[0]+i/((16-i)*x[1]+Math.min(i,16-i)*x[2]))); return ret; }; x0 = numeric.uncmin(f,[1,1,1]).solution
[0.08241,1.133,2.344]
> y = [0.0009,0.0044,0.0175,0.0540,0.1295,0.2420,0.3521,0.3989,0.3521,0.2420,0.1295,0.0540,0.0175,0.0044,0.0009]; f = function(x) { var ret = 0,i; for(i=1;i<=15;i++) ret+=sqr(x[0]*Math.exp(-x[1]*sqr((8-i)/2-x[2])/2)-y[i-1]); return ret; }; x0 = numeric.div(numeric.round(numeric.mul(numeric.uncmin(f,[1,1,1]).solution,1000)),1000)
[0.399,1,0]
> f = function(x) { return sqr(x[0]+10*x[1])+5*sqr(x[2]-x[3])+sqr(sqr(x[1]-2*x[2]))+10*sqr(x[0]-x[3]); }; x0 = numeric.div(numeric.round(numeric.mul(numeric.uncmin(f,[3,-1,0,1]).solution,1e5)),1e5)
[0,0,0,0]
> f = function(x) { return sqr(10*(x[1]-x[0]*x[0]))+sqr(1-x[0])+90*sqr(x[3]-x[2]*x[2])+sqr(1-x[2])+10*sqr(x[1]+x[3]-2)+0.1*sqr(x[1]-x[3]); }; x0 = numeric.uncmin(f,[-3,-1,-3,-1]).solution
[1,1,1,1]
> y = [0.1957,0.1947,0.1735,0.1600,0.0844,0.0627,0.0456,0.0342,0.0323,0.0235,0.0246]; u = [4,2,1,0.5,0.25,0.167,0.125,0.1,0.0833,0.0714,0.0625]; f = function(x) { var ret=0, i; for(i=0;i<11;++i) ret += sqr(y[i]-x[0]*(u[i]*u[i]+u[i]*x[1])/(u[i]*u[i]+u[i]*x[2]+x[3])); return ret; }; x0 = numeric.uncmin(f,[0.25,0.39,0.415,0.39]).solution
[ 0.1928, 0.1913, 0.1231, 0.1361]
> y = [0.844,0.908,0.932,0.936,0.925,0.908,0.881,0.850,0.818,0.784,0.751,0.718,0.685,0.658,0.628,0.603,0.580,0.558,0.538,0.522,0.506,0.490,0.478,0.467,0.457,0.448,0.438,0.431,0.424,0.420,0.414,0.411,0.406]; f = function(x) { var ret=0, i; for(i=0;i<33;++i) ret += sqr(y[i]-(x[0]+x[1]*Math.exp(-10*i*x[3])+x[2]*Math.exp(-10*i*x[4]))); return ret; }; x0 = numeric.uncmin(f,[0.5,1.5,-1,0.01,0.02]).solution
[ 0.3754, 1.936, -1.465, 0.01287, 0.02212]
> f = function(x) { var ret=0, i,ti,yi,exp=Math.exp; for(i=1;i<=13;++i) { ti = 0.1*i; yi = exp(-ti)-5*exp(-10*ti)+3*exp(-4*ti); ret += sqr(x[2]*exp(-ti*x[0])-x[3]*exp(-ti*x[1])+x[5]*exp(-ti*x[4])-yi); } return ret; }; x0 = numeric.uncmin(f,[1,2,1,1,1,1],1e-14).solution; f(x0)<1e-20;
true
There are optional parameters to numeric.uncmin(f,x0,tol,gradient,maxit,callback). The iteration stops when
the gradient or step size is less than the optional parameter tol. The gradient() parameter is a function that computes the
gradient of f(). If it is not provided, a numerical gradient is used. The iteration stops when
maxit iterations have been performed. The optional callback() parameter, if provided, is called at each step:
> z = []; cb = function(i,x,f,g,H) { z.push({i:i, x:x, f:f, g:g, H:H }); }; x0 = numeric.uncmin(function(x) { return Math.cos(2*x[0]); },[1],1e-10,function(x) { return [-2*Math.sin(2*x[0])]; },100,cb);
{solution: [1.571],
f: -1,
gradient: [2.242e-11],
invHessian: [[0.25]],
iterations: 6,
message: "Newton step smaller than tol"}
> z
[{i:0, x:[1 ], f:-0.4161, g: [-1.819 ] , H:[[1 ]] },
{i:2, x:[1.909], f:-0.7795, g: [ 1.253 ] , H:[[0.296 ]] },
{i:3, x:[1.538], f:-0.9979, g: [-0.1296 ] , H:[[0.2683]] },
{i:4, x:[1.573], f:-1 , g: [ 9.392e-3] , H:[[0.2502]] },
{i:5, x:[1.571], f:-1 , g: [-6.105e-6] , H:[[0.25 ]] },
{i:6, x:[1.571], f:-1 , g: [ 2.242e-11] , H:[[0.25 ]] }]
Solving ODEs
The function numeric.dopri() is an implementation of the Dormand-Prince-Runge-Kutta integrator with
adaptive time-stepping:
> sol = numeric.dopri(0,1,1,function(t,y) { return y; })
{ x: [ 0, 0.1, 0.1522, 0.361, 0.5792, 0.7843, 0.9813, 1],
y: [ 1, 1.105, 1.164, 1.435, 1.785, 2.191, 2.668, 2.718],
f: [ 1, 1.105, 1.164, 1.435, 1.785, 2.191, 2.668, 2.718],
ymid: [ 1.051, 1.134, 1.293, 1.6, 1.977, 2.418, 2.693],
iterations:8,
events:,
message:""}
The return value sol contains the x and y values of the solution.
If you need to know the value of the solution between the given x values, use the function
sol.at(), which uses the extra information contained in sol.ymid and sol.f to
produce approximations between these points:
> sol.at([0.3,0.7])
[1.35,2.014]
The integrator is also able to handle vector equations. This is a harmonic oscillator:
> numeric.dopri(0,10,[3,0],function (x,y) { return [y[1],-y[0]]; }).at([0,0.5*Math.PI,Math.PI,1.5*Math.PI,2*Math.PI])
[[ 3, 0],
[ -9.534e-8, -3],
[ -3, 2.768e-7],
[ 3.63e-7, 3],
[ 3, -3.065e-7]]
Van der Pol:
> numeric.dopri(0,20,[2,0],function(t,y) { return [y[1], (1-y[0]*y[0])*y[1]-y[0]]; }).at([18,19,20])
[[ -1.208, 0.9916],
[ 0.4258, 2.535],
[ 2.008, -0.04251]]
You can also specify a tolerance, a maximum number of iterations and an event function. The integration stops if
the event function goes from negative to positive.
> sol = numeric.dopri(0,2,1,function (x,y) { return y; },1e-8,100,function (x,y) { return y-1.3; });
{ x: [ 0, 0.0181, 0.09051, 0.1822, 0.2624],
y: [ 1, 1.018, 1.095, 1.2, 1.3],
f: [ 1, 1.018, 1.095, 1.2, 1.3],
ymid: [ 1.009, 1.056, 1.146, 1.249],
iterations: 5,
events: true,
message: "" }
Note that at x=0.1822, the event function value was -0.1001 while at x=0.2673, the event
value was 6.433e-3. The integrator thus terminated at x=0.2673 instead of continuing until the end
of the integration interval.
Events can also be vector-valued:
> sol = numeric.dopri(0,2,1,function(x,y) { return y; },undefined,50,function(x,y) { return [y-1.5,Math.sin(y-1.5)]; });
{ x: [ 0, 0.2, 0.4055],
y: [ 1, 1.221, 1.5],
f: [ 1, 1.221, 1.5],
ymid: [1.105, 1.354],
iterations: 2,
events: [true,true],
message: ""}
Seedrandom (David Bau)
The object numeric.seedrandom is based on
David Bau's seedrandom.js.
This small library can be used to create better pseudorandom numbers than Math.random() which can
furthermore be "seeded".
> numeric.seedrandom.seedrandom(3); numeric.seedrandom.random()
0.7569
> numeric.seedrandom.random()
0.6139
> numeric.seedrandom.seedrandom(3); numeric.seedrandom.random()
0.7569
For performance reasons, numeric.random() uses the default Math.random(). If you want to use
the seedrandom version, just do Math.random = numeric.seedrandom.random. Note that this may slightly
decrease the performance of all Math operations.