Numerics

RepLAB provides the following multi-dimensional array types as an extension to the types provided by MATLAB/Octave.

  • cyclotomic implements exact computation on the cyclotomic field. This enables RepLAB to handle any irreducible representation of a finite group, if care is taken to express it over the cyclotomics (which is always possible).

  • H implements matrices of quaternions using a pair of complex matrices. Quaternion matrices are used to decompose real representations of quaternion-type.

cyclotomic

class replab.cyclotomic(array, size_)

Multidimensional array with elements in the cyclotomic field

Requires Java and the cyclolab support library.

The implementation is pretty slow, but values are exact.

We can construct cyclotomic arrays in various ways.

They can be constructed from double floating-point arrays. Note that only fractions with a power-of-two denominator can be represented exactly.

Example

>>> replab.cyclotomic(1/2) % doctest: +cyclotomic
    1/2
>>> replab.cyclotomic(1/3)
    6004799503160661/18014398509481984

Note that passing a cyclotomic array as a parameter is the identity operation. This is useful when converting number types.

Example

>>> c = replab.cyclotomic(1) % doctest: +cyclotomic
    1
>>> replab.cyclotomic(c)
    1

Cyclotomics can also be constructed from vpis:

Example

>>> c = replab.cyclotomic(vpi('100000000000000')) % doctest: +cyclotomic
    100000000000000

They can also be constructed from strings. A vector/matrix syntax is also available, but row coefficients must always be separated by commas.

Example

>>> replab.cyclotomic('2/3') % doctest: +cyclotomic
    2/3
>>> replab.cyclotomic('[1, 0; 0, 1]')
    1  0
    0  1

The string syntax accepts:

  • integers,

  • the four operations +, -, *, /,

  • powers ^ with integer exponents,

  • parentheses,

  • roots of unity as in Gap System (E(n) is exp(2*i*pi/n)),

  • square roots of rational numbers,

  • sines and cosines of rational multiples of pi.

Example

>>> replab.cyclotomic('sqrt(2)') % doctest: +cyclotomic
    E(8)-E(8)^3
>>> replab.cyclotomic('cos(-pi/3)') % doctest: +cyclotomic
    -1/2
>>> replab.cyclotomic('(3 + 2/3)^3') % doctest: +cyclotomic
    1331/27
>>> replab.cyclotomic('E(3)^2') % doctest: +cyclotomic
    E(3)^2

Finally, the constructor also accepts heterogenous cell arrays, where the types above can be mixed and matched (except that strings must represent scalars, not matrices/vectors).

Example

>>> replab.cyclotomic({'1' vpi(2); '1/2' 1}) % doctest: +cyclotomic
      1    2
     1/2   1

A variety of static methods is available.

Example

>>> I = replab.cyclotomic.eye(3) % doctest: +cyclotomic
    1  0  0
    0  1  0
    0  0  1

Own members

data_

1D array of cyclotomic numbers, elements are stored as in array(:)

Type

cyclo.Cyclo[]

size_

Shape

Type

integer(1,*)

static E(orders)

Returns roots of unity

Parameters

n (integer(*,*)) – Root orders

Returns

The value exp(2i*pi/n)

Return type

cyclotomic

static approximate(lowerBounds, upperBounds)

Constructs the simplest rational approximation within an interval, for coefficients of a matrix

Example

>>> replab.cyclotomic.approximate(3.141592, 3.141593) % doctest: +cyclotomic
  355/113
Parameters
  • lowerBounds (double(...)) – Interval lower bounds

  • upperBounds (double(...)) – Interval upper bounds

Returns

The simplest rational approximation, coefficient by coefficient

Return type

cyclotomic

blkdiag(lhs, varargin)

Block diagonal concatenation of matrix input arguments

Compatible with MATLAB’s blkdiag

conj(self)

Complex conjugation

Fully compatible with MATLAB’s conj

ctranspose(self)

Complex conjugate transpose

Fully compatible with MATLAB’s ctranspose

cyclo(self, varargin)

Returns the cyclo.Cyclo object at the given index in this multidimensional array

Arguments are integer indices, as passed to subsref

cyclotomic(array, size_)

Constructs a cyclotomic array from an array of coefficients

Parameters

array – Coefficients

data(self)

Returns the java array containing the flat data of this array

diag(self, varargin)

Diagonal matrices and diagonals of a matrix

Fully compatible with MATLAB’s diag

disp(self)

Standard display method

dot(lhs, rhs)

Vector dot product (incomplete)

Only supports the case where both lhs and rhs are vectors, and does not accept a dim argument.

double(self)

Conversion to floating-point double

doubleApproximation(self)

Conversion to floating-point double with estimation of the error

eq(lhs, rhs)

Equality test, returns a multidimensional array

static eye(n)

Constructs the n x n identity matrix (incomplete)

Does not support the syntax with two arguments.

Parameters

n (integer) – Matrix size

Returns

The constructed matrix

Return type

cyclotomic

find(self, varargin)

Find indices of nonzero elements

Compatible with MATLAB’s find

static fromRationals(numerators, denominators)

Constructs a cyclotomic matrix of rational numbers

Example

>>> replab.cyclotomic.fromRationals([1 1; 1 1], [2 3; 3 2]) % doctest: +cyclotomic
    1/2  1/3
    1/3  1/2
Parameters
  • numerators (double(*,*)) – Coefficient numerators

  • denominators (double(*,*)) – Coefficient denominators

Returns

The constructed matrix

Return type

cyclotomic

full(self)

Placeholder as cyclotomic arrays are always dense

galois(self, ord)

Action of the Galois group

hash(self)

Returns a multidimensional array of hash codes

Example

>>> hash(replab.cyclotomic.zeros(1)) == -1579025880 % doctest: +cyclotomic
    1
Returns

Integer array of hash codes

Return type

integer(…)

horzcat(lhs, varargin)

Horizontal concatenation (incomplete)

Incomplete implementation: only works on vectors and matrices

imag(self)

Returns the imaginary part

inv(self)

Matrix inverse

Requires that the matrix is full rank

ipermute(self, order)

Inverse permute array dimensions

Fully compatible with MATLAB’s permute

isequal(lhs, rhs)

Equality test, returns a scalar

isequaln(lhs, rhs)

Equality test, returns a scalar

Same as isequal as we do not have NaN

isrational(self)

Returns which coefficients are rational

iswhole(self)

Returns which coefficients are integers

length(self)

Length of vector

Equivalent to max(size(X)) if isempty(X) is false, otherwise 0.

Example

>>> length(replab.cyclotomic.zeros(10)) % doctest: +cyclotomic
    10
Returns

Vector length

Return type

integer

lu(self)

Returns the LU decomposition of this matrix

The m x n matrix must have m >= n.

It returns triangular matrices L, U and a permutation vector p such that L*U == self(:,p)

Returns

  • L (cyclotomic (*,*)) – Lower triangular matrix

  • U (cyclotomic (*,*)) – Upper triangular matrix

  • p (integer(1,*)) – Permutation vector

minus(lhs, rhs)

Standard - operator

mpower(self, e)

Matrix power

mrdivide(lhs, rhs)

Division

Only supports the case where the right hand side is scalar

mtimes(lhs, rhs)

Matrix multiplication

Support both m x n by n x p matrix multiplication, and the case where one of the arguments is a scalar

ndims(self)

Number of dimensions

Returns

Number of dimensions

Return type

integer

ne(self, rhs)

(Non-)equality test

null(self)

Computes the null space of a cyclotomic matrix

Example

>>> M = replab.cyclotomic([1 3 0; -2 -6 0; 3 9 6]); % doctest: +cyclotomic
>>> null(M)'
    -3  1  0
Returns

The matrix null space

Return type

cyclotomic

num2str(self)

Conversion to string (incomplete)

Only supports matrices or vectors.

Does not support format specifiers or a specified precision.

Returns

Possibly multiline string representation of the matrix

Return type

charstring

numel(self)

Number of elements

Returns

Number of elements

Return type

integer

static ones(varargin)

Constructs a cyclotomic array filled with ones

Returns

The constructed matrix

Return type

cyclotomic

permute(self, order)

Permute array dimensions

Fully compatible with MATLAB’s permute

plus(lhs, rhs)

Standard + operator

prod(self, dim)

Product of elements (incomplete)

Does not support arrays with more than two dimensions, and does not support extra options (type, nanflag).

static rand(varargin)

Random cyclotomic matrix

Parameters
  • maximalNumberOfTerms – Maximal number of terms

  • field – Whether to sample from real or complex cyclotomics

  • maximalOrder – Maximal root-of-unity order

  • maximalDenominator – Largest denominator for rational numbers

  • maximalAbsoluteValue – Approximate largest absolute value

Kwtype maximalNumberOfTerms

integer

Kwtype field

‘R’, ‘C’, ‘Q’

Kwtype maximalOrder

integer

Kwtype maximalDenominator

integer

Kwtype maximalAbsoluteValue

integer

rdivide(lhs0, rhs0)

Pointwise / operator

real(self)

Returns the real part

reshape(self, varargin)

Reshape array

Fully compatible with MATLAB’s reshape

rref(self)

Computes the row reduced echelon form

Returns

The row reduced echelon form

Return type

cyclotomic

size(self, dim)

Size of array (incomplete)

Incomplete implementation: only supports the forms s = size(array) and s = size(array, i) where i is an integer scalar.

Example

>>> size(replab.cyclotomic.eye(2)) % doctest: +cyclotomic
    2     2
Parameters

dim (integer, optional) – Dimension(s) to get the size of

Returns

Size

Return type

integer or integer(1,*)

static sparse(I, J, K, nR, nC)

Constructs a 2D cyclotomic array using sparse matrix syntax (incomplete)

Only supports the five argument syntax

Returns

The constructed matrix

Return type

cyclotomic

sqrt(self)

Square root

Requires that the argument contains rational coefficients, otherwise compatible with MATLAB’s sqrt

static sqrtRational(num, den)

Returns the square root of a rational number

Example

>>> replab.cyclotomic.sqrtRational(2) % doctest: +cyclotomic
    E(8)-E(8)^3
>>> replab.cyclotomic.sqrtRational(1,2)
    E(8)/2-E(8)^3/2
Parameters
  • num (integer) – Numerator

  • den (integer, optional) – Denominator, default value 1

Returns

The value sqrt(num/den)

Return type

cyclotomic

subsasgn(self, s, val)

Assignment

Compatible with MATLAB’s subsasgn

subsref(self, s)

Indexing

Compatible with MATLAB’s subsref

sum(self, dim)

Sum of elements (incomplete)

Does not support arrays with more than two dimensions, and does not support extra options (type, nanflag).

times(lhs0, rhs0)

Pointwise * operator

toCellOfStrings(self)

Returns a multi-dimensional cell array of strings corresponding to this cyclotomic array

trace(self)

Sum of diagonal elements

Fully compatible with MATLAB’s trace

transpose(self)

Transpose

Fully compatible with MATLAB’s transpose

uminus(self)

Unary minus

vertcat(lhs, varargin)

Vertical concatenation

static zeros(varargin)

Constructs a cyclotomic array filled with zeros

Returns

The constructed matrix

Return type

cyclotomic

H

class replab.H(M1, Mi, Mj, Mk)

Quaternion matrix type

A quaternion matrix can be expressed using four real matrices A, B, C and D:

\(Q = A + B i + C j + D k\) where \(i^2 = j^2 = k^2 = i j k = -1\) and \(i j = k\), \(j k = i\) and \(k i = j\).

We choose to encode our quaternions using two matrices \(X = A + B i\) and \(Y = C + D i\), where we now identify \(i\) with the complex imaginary unit.

Thus \(Q = X + Y j\).

The class supports sparse components for efficiency.

Own members

X

First part

Y

Second part

H(M1, Mi, Mj, Mk)

Quaternion constructor

Constructs a quaternion as: `` Q = M1 + Mi * i + Mj * j + Mk * k `` where i, j, k are the quaternion units.

Some matrices can be complex, in which case the complex imaginary unit is identified with i.

Parameters
  • M1 (double(d1,d2)) – Must be provided

  • Mi (double(d1,d2) or []) – Matrix to be multiplied by i

  • Mj (double(d1,d2) or []) – Matrix to be multiplied by j

  • Mk (double(d1,d2) or []) – Matrix to be multiplied by k

abs(self)

Computes the absolute value, array coefficient by coefficient

Returns

Magnitude of each quaternion coefficient

Return type

double

abs2(self)

Computes the square of the absolute value, array coefficient by coefficient

Returns

Square of the magnitude of each quaternion coefficient

Return type

double

static decode(Z)

Decodes the given complex matrix into a quaternion matrix

Parameters

Z (double(*,*)) – Complex matrix encoding a quaternion algebra, see encode

Returns

Quaternion matrix

Return type

replab.H

static decompose(Q)

Decomposes the first and second part from the given matrix, interpreted as a quaternion matrix

The result of this call is equivalent to tmp = replab.H(Q); X = tmp.X; Y = tmp.Y, but it avoids unnecessary construction of temporary objects.

Parameters

Q (replab.H or double multidimensional array) – Quaternion to decompose

Returns

  • X (double multidimensional array) – First part of the quaternion matrix

  • Y (double multidimensional array) – Second part of the quaternion matrix

disp(self)

Standard display method

static encode(Q)

Encodes the given matrix, interpreted as a quaternion, into a complex matrix

Uses the encoding of a quaternion q = a + b i + c j + d k as [a + i*b, c + i*d; -c + i*d, a - i*b]

Parameters

Q (H) – Quaternion matrix

Returns

Complex matrix encoding the given quaternion matrix

Return type

double(*,*)

eps(self)

Computes the spacing of floating point numbers

static j()

Quaternion unit j

Returns

Quaternion scalar matrix

Return type

H

static k()

Quaternion unit k

Returns

Quaternion scalar matrix

Return type

H

mldivide(Q1, Q2)

Backslash or left matrix divide

mrdivide(Q1, Q2)

Slash or right matrix divide

part1(self)

Returns the real part of the quaternion

Returns

Real part

Return type

double(…)

parti(self)

Returns the i part of the quaternion

Returns

Part

Return type

double(…)

partj(self)

Returns the j part of the quaternion

Returns

Part

Return type

double(…)

partk(self)

Returns the k part of the quaternion

Returns

Part

Return type

double(…)

static randn(varargin)

Returns a normally distributed quaternion multidimensional array

Parameters

varargin – Arguments passed to randn

Returns

Quaternion multidimensional array

Return type

H

sqrt(self)

Quaternion element-by-element square root

See https://math.stackexchange.com/questions/382431/square-roots-of-quaternions

toStringCell(self)

Computes a string representation of the elements of the array

Returns

String representation

Return type

cell(…) of charstring