Skip to content

Breeze Linear Algebra

Tim Bezemer edited this page Sep 25, 2020 · 16 revisions

Breeze Linear Algebra

Core Concepts

Compared to other numerical computing environments, Breeze matrices default to column major ordering, like Matlab, but indexing is 0-based, like Numpy. Breeze has as its core concepts matrices and column vectors. Row vectors are normally stored as matrices with a single row. This allows for greater type safety with the downside that conversion of row vectors to column vectors is performed using a transpose-slice (a.t(::,0)) instead of a simple transpose (a.t).

UFuncs are very important in Breeze. Once you get a feel for the syntax (i.e. what's in this section), it might be worthwhile to read the first half of the UFunc page. (You can skip the last half that involves implementing your own UFuncs...until you're ready to contribute to Breeze!)

Quick Reference

The following table assumes that Numpy is used with from numpy import * and Breeze with

import breeze.linalg._
import breeze.numerics._

Creation

Operation Breeze Matlab Numpy
Zeroed matrix val c = DenseMatrix.zeros[Double](n,m) c = zeros(n,m) c = zeros((n,m))
Zeroed vector val a = DenseVector.zeros[Double](n) a = zeros(n) a = zeros(n)
Vector of ones val a = DenseVector.ones[Double](n) a = ones(n) a = ones(n)
Vector of particular number val b = DenseVector.fill(n){5.0} a = ones(n) * 5 a = ones(n) * 5
Identity matrix DenseMatrix.eye[Double](n) eye(n) eye(n)
Diagonal matrix diag(DenseVector(1.0,2.0,3.0)) diag([1 2 3]) diag((1,2,3))
Matrix inline creation val a = DenseMatrix((1.0,2.0), (3.0,4.0)) a = [1 2; 3 4] a = array([ [1,2], [3,4] ])
Column vector inline creation val a = DenseVector(1,2,3,4) a = [1 2 3 4] a=array([1,2,3,4])
Row vector inline creation val a = DenseVector(1,2,3,4).t a = [1,2,3,4]' a = array([1,2,3]).reshape(-1,1)

Indexing and Slicing

Operation Breeze Matlab Numpy
Basic Indexing a(0,1) a(1,2) a[0,1]
Extract subset of vector a(1 to 4) or a(1 until 5) or a.slice(1,4) a(2:5) a[1:4]
(negative steps) a(5 to 0 by -1) a(6:-1:1) a[5:0:-1]
(tail) a(1 to -1) a(2:end) a[1:]
(last element) a( -1 ) a(end) a[-1]
Extract column of matrix a(::, 2) a(:,3) a[:,2]

Other Manipulation

Operation Breeze Matlab Numpy
Reshaping a.reshape(3, 2) reshape(a, 3, 2) a.reshape(3,2)
Flatten matrix a.toDenseVector (Makes copy) a(:) a.flatten()
Copy lower triangle lowerTriangular(a) tril(a) tril(a)
Copy upper triangle upperTriangular(a) triu(a) triu(a)
Create view of matrix diagonal diag(a) NA diagonal(a) (Numpy >= 1.9)
Vector Assignment to subset a(1 to 4) := 5.0 a(2:5) = 5 a[1:4] = 5
Vector Assignment to subset a(1 to 4) := DenseVector(1.0,2.0,3.0) a(2:5) = [1 2 3] a[1:4] = array([1,2,3])
Matrix Assignment to subset a(1 to 3,1 to 3) := 5.0 a(2:4,2:4) = 5 a[1:3,1:3] = 5
Matrix Assignment to column a(::, 2) := 5.0 a(:,3) = 5 a[:,2] = 5
Matrix vertical concatenate DenseMatrix.vertcat(a,b) [a ; b] vstack((a,b))
Matrix horizontal concatenate DenseMatrix.horzcat(d,e) [a , b] hstack((a,b))
Vector concatenate DenseVector.vertcat(a,b) [a b] concatenate((a,b))

Operations

Operation Breeze Matlab Numpy
Elementwise addition a + b a + b a + b
Elementwise multiplication a :* b a .* b a * b
Elementwise exponentiation a ^:^ b a .^ b power(a, b)
Elementwise comparison a :< b a < b (gives matrix of 1/0 instead of true/false) a < b
Elementwise equals a :== b a == b (gives matrix of 1/0 instead of true/false) a == b
Inplace addition a :+= 1.0 a += 1 a += 1
Inplace elementwise multiplication a :*= 2.0 a *= 2 a *= 2
Vector dot product a dot b, a.t * b dot(a,b) dot(a,b)
Elementwise max max(a) max(a) a.max()
Elementwise argmax argmax(a) [v i] = max(a); i a.argmax()

†: new in 0.7-SNAPSHOT

Sum

Operation Breeze Matlab Numpy
Elementwise sum sum(a) sum(sum(a)) a.sum()
Sum down each column (giving a row vector) sum(a, Axis._0) or sum(a(::, *)) sum(a) sum(a,0)
Sum across each row (giving a column vector) sum(a, Axis._1) or sum(a(*, ::)) sum(a') sum(a,1)
Trace (sum of diagonal elements) trace(a) trace(a) a.trace()
Cumulative sum accumulate(a) cumsum(a) a.cumsum()

Boolean Operators

Operation Breeze Matlab Numpy
Elementwise and a :& b a && b a and b
Elementwise or `a : b` `a
Elementwise not !a ~a not a
True if any element is nonzero any(a) any(a) n/a
True if all elements are nonzero all(a) all(a) n/a

Linear Algebra Functions

Operation Breeze Matlab Numpy
Linear solve a \ b a \ b linalg.solve(a,b)
Transpose a.t a' a.conj.transpose()
Determinant det(a) det(a) linalg.det(a)
Inverse inv(a) inv(a) linalg.inv(a)
Pseudoinverse pinv(a) pinv(a) linalg.pinv(a)
Norm norm(a) norm(a) norm(a)
Eigenvalues (Symmetric) eigSym(a) [v,l] = eig(a) linalg.eig(a)[0]
Eigenvalues val (er, ei, _) = eig(a) (Separate real & imaginary part) eig(a) linalg.eig(a)[0]
Eigenvectors eig(a)._3 [v,l] = eig(a) linalg.eig(a)[1]
Singular values val (u,s,v) = svd(a) svd(a) linalg.svd(a)
Rank rank(a) rank(a) rank(a)
Vector length a.length size(a) a.size
Matrix rows a.rows size(a,1) a.shape[0]
Matrix columns a.cols size(a,2) a.shape[1]

Rounding and Signs

Operation Breeze Matlab Numpy
Round round(a) round(a) around(a)
Ceiling ceil(a) ceil(a) ceil(a)
Floor floor(a) floor(a) floor(a)
Sign signum(a) sign(a) sign(a)
Absolute Value abs(a) abs(a) abs(a)

Constants

Operation Breeze Matlab Numpy
Not a Number NaN or nan NaN nan
Infinity Inf or inf Inf inf
Pi Constants.Pi pi math.pi
e Constants.E exp(1) math.e

Complex numbers

If you make use of complex numbers, you will want to include a breeze.math._ import. This declares a i variable, and provides implicit conversions from Scala’s basic types to complex types.

Operation Breeze Matlab Numpy
Imaginary unit i i z = 1j
Complex numbers 3 + 4 * i or Complex(3,4) 3 + 4i z = 3 + 4j
Absolute Value abs(z) or z.abs abs(z) abs(z)
Real Component z.real real(z) z.real
Imaginary Component z.conjugate conj(z) z.conj() or z.conjugate()

Numeric functions

Breeze contains a fairly comprehensive set of special functions under the Breeze.numerics._ import. These functions can be applied to single elements, vectors or matrices of Doubles. This includes versions of the special functions from scala.math that can be applied to vectors and matrices. Any function acting on a basic numeric type can “vectorized”, to a UFunc function, which can act elementwise on vectors and matrices:

val v = DenseVector(1.0,2.0,3.0)
exp(v) // == DenseVector(2.7182818284590455, 7.38905609893065, 20.085536923187668)

UFuncs can also be used in-place on Vectors and Matrices:

val v = DenseVector(1.0,2.0,3.0)
exp.inPlace(v) // == DenseVector(2.7182818284590455, 7.38905609893065, 20.085536923187668)

See UFunc for more information.

Here is a (non-exhaustive) list of UFuncs in Breeze:

Trigonometry

  • sin, sinh, asin, asinh
  • cos, cosh, acos, acosh
  • tan, tanh, atan, atanh
  • atan2
  • sinc(x) == sin(x)/x
  • sincpi(x) == sinc(x * Pi)

Logarithm, Roots, and Exponentials

  • log, exp log10
  • log1p, expm1
  • sqrt, sbrt
  • pow

Gamma Function and its cousins

The gamma function is the extension of the factorial function to the reals. Numpy needs from scipy.special import * for this and subsequent sections.

Operation Breeze Matlab Numpy
Gamma function exp(lgamma(a)) gamma(a) gamma(a)
log Gamma function lgamma(a) gammaln(a) gammaln(a)
Incomplete gamma function gammp(a, x) gammainc(a, x) gammainc(a, x)
Upper incomplete gamma function gammq(a, x) gammainc(a, x, tail) gammaincc(a, x)
derivative of lgamma digamma(a) psi(a) polygamma(0, a)
derivative of digamma trigamma(a) psi(1, a) polygamma(1, a)
nth derivative of digamma na psi(n, a) polygamma(n, a)
Log Beta function lbeta(a,b) betaln(a, b) betaln(a,b)
Generalized Log Beta function lbeta(a) na na

Error Function

The error function...

Operation Breeze Matlab Numpy
error function erf(a) erf(a) erf(a)
1 - erf(a) erfc(a) erfc(a) erfc(a)
inverse error function erfinv(a) erfinv(a) erfinv(a)
inverse erfc erfcinv(a) erfcinv(a) erfcinv(a)

Other functions

Operation Breeze Matlab Numpy
logistic sigmoid sigmoid(a) na expit(a)
Indicator function I(a) not needed where(cond, 1, 0)

Signal Analysis

Signal analysis functions, including correlate and convolve live in scala.breeze, and are documented in Breeze Signal Processing

Broadcasting

Sometimes we want to apply an operation to every row or column of a matrix, as a unit. For instance, you might want to compute the mean of each row, or add a vector to every column. Adapting a matrix so that operations can be applied columnwise or rowwise is called broadcasting. Languages like R and numpy automatically and implicitly do broadcasting, meaning they won’t stop you if you accidentally add a matrix and a vector. In Breeze, you have to signal your intent using the broadcasting operator *. The * is meant to evoke “foreach” visually. Here are some examples:

    val dm = DenseMatrix((1.0,2.0,3.0),
                         (4.0,5.0,6.0))

    val res = dm(::, *) + DenseVector(3.0, 4.0)
    assert(res === DenseMatrix((4.0, 5.0, 6.0), (8.0, 9.0, 10.0)))

    res(::, *) := DenseVector(3.0, 4.0)
    assert(res === DenseMatrix((3.0, 3.0, 3.0), (4.0, 4.0, 4.0)))

    val m = DenseMatrix((1.0, 3.0), (4.0, 4.0))
    assert(mean(m(*, ::)) === DenseVector(2.0, 4.0))
    assert(mean(m(::, *)) === DenseMatrix((2.5, 3.5)))

Elementwise special functions

Breeze contains a fairly comprehensive set of special functions under the Breeze.numerics._ import. These functions can be applied to single elements, vectors or matrices of Doubles. This includes versions of the special functions from scala.math that can be applied to vectors and matrices. Any function acting on a basic numeric type can “vectorized”, to a UFunc function, which can act elementwise on vectors and matrices:

val v = DenseVector(1.0,2.0,3.0)
exp(v) // == DenseVector(2.7182818284590455, 7.38905609893065, 20.085536923187668)

UFuncs can also be used in-place on Vectors and Matrices:

val v = DenseVector(1.0,2.0,3.0)
exp(v) // == DenseVector(2.7182818284590455, 7.38905609893065, 20.085536923187668)

Some UFuncs are more like “reduction” functions than like “vectorized functions.” For instance, sum computes the sum of the entries in a matrix or vector.

val r = DenseMatrix.rand(5,5)

// sum all elements
sum(r):Double 

sum and mean can also operate on the rows or columns of a matrix independently, using broadcasting:

// mean of each row into a single column
mean(r(*, ::))

// sum of each column into a single row
sum(r(::, *))

The UFunc trait is similar to numpy’s concept of ufuncs. See UFunc for more information on Breeze UFuncs.

Casting and type safety

Compared to Numpy and Matlab, Breeze requires you to be more explicit about the types of your variables. When you create a new vector for example, you must specify a type (such as in DenseVector.zeros[Double](n)) in cases where a type can not be inferred automatically. Automatic inference will occur when you create a vector by passing its initial values in (DenseVector). A common mistake is using integers for initialisation (i.e. DenseVector), which would give a matrix of integers instead of doubles. Both Numpy and Matlab would default to doubles instead.

Breeze will not convert integers to doubles for you in most expressions. Simple operations like a :+ 3 when a is a DenseVector[Double] will not compile. Breeze provides a convert function, which can be used to explicitly cast. You can also use v.mapValues(_.toDouble).

Casting

Operation Breeze Matlab Numpy
Convert to Int convert(a, Int) int(a) a.astype(int)

Performance

Breeze uses netlib-java for its core linear algebra routines. This includes all the cubic time operations, matrix-matrix and matrix-vector multiplication. Special efforts are taken to ensure that arrays are not copied.

Netlib-java will attempt to load system optimised BLAS/LAPACK if they are installed, falling back to the reference natives, falling back to pure Java. Set your logger settings to ALL for the com.github.fommil.netlib package to check the status, and to com.github.fommil.jniloader for a more detailed breakdown. Read the netlib-java project page for more details.

Currently vectors and matrices over types other than Double, Float and Int are boxed, and will typically be a lot slower. If you find yourself needing other AnyVal types like Long or Short, please ask on the list about possibly adding support for them.