Equivariant spaces

Equivariant spaces describe linear map between representations that preserve the group structure.

The hierarchy of equivariant spaces closely mirrors the hierarchy of representations Rep.

  • Equivariant describes the vector space of equivariant linear maps between two representations of the same group.

  • SubEquivariant describes a subspace of an equivariant space induced by subspaces of source/target representations.

  • IsotypicEquivariant describes a subspace of an equivariant space corresponding to a pair of isotypic subspaces of the source/target representations

  • IrreducibleEquivariant is the most interesting class, as it describes the equivariant space in its block form, taking in account the decomposition of the source and target representation into irreps.

Equivariant

class replab.Equivariant(repR, repC, special)

Bases: replab.Domain

Describes a vector space of group-equivariant matrices

Let repC and repR be two representations of the same group G.

This describes the set of matrices X such that repR.image(g) * X = X * repC.image(g)

See Proposition 4 of J.-P. Serre, Linear Representations of Finite Groups (Springer, 1977).

There are several special cases of equivariant spaces which RepLAB uses extensively; these are parameterized by a single representation rho from which the equivariant space is constructed. They are distinguished by the special property.

Special name

repR

repC

Additional constraint

antilinear

rho

conj(rho)

bilinear

dual(rho)

rho

commutant

rho

rho

hermitian

conj(dual(rho))

rho

X = X’

sesquilinear

conj(dual(rho))

rho

symmetric

dual(rho)

rho

X = X.’

trivialRows

trivial (d = rho.dimension)

rho

trivialCols

rho

trivial (d = rho.dimension)

When rho is unitary:

  • commutant and sesquilinear spaces are identical,

  • antilinear and bilinear space are identical.

When rho is real:

  • antilinear and commutant spaces are identical,

  • bilinear and sesquilinear spaces are identical,

  • hermitian and symmetric spaces are identical.

Own members

domain

Domain, real or complex matrices

Type

replab.Domain

field

Field of the vector space real (R) or complex x(C)

Type

{‘R’, ‘C’}

group

Group being represented

Type

replab.CompactGroup

nC

Column size

Type

integer

nR

Row size

Type

integer

repC

Representation of column space

Type

replab.Rep

repR

Representation of row space

Type

replab.Rep

special

Whether the equivariant space has special structure, see comment in replab.Equivariant

Type

charstring

decomposition(self)

Returns the decomposition of this equivariant space

Returns

The equivariant space

Return type

IrreducibleEquivariant

isExact(self)

Returns whether this equivariant space can compute exact projection

Returns

True if the call self.projection(X, 'exact') always succeeds

Return type

logical

static make(repR, repC, varargin)

Returns the space of equivariant linear maps between two representations

The equivariant vector space contains the matrices X such that

repC.image(g) * X = X * repR.image(g)

Parameters
  • repR (Rep) – Representation on the target/row space

  • repC (Rep) – Representation on the source/column space

  • special – Special structure if applicable, see Equivariant, default: ‘’

  • type – Whether to obtain an exact equivariant space, default ‘double’ (‘double’ and ‘double/sparse’ are equivalent)

Kwtype special

charstring, optional

Kwtype type

‘exact’, ‘double’ or ‘double/sparse’, optional

Returns

The equivariant vector space

Return type

Equivariant

project(self, X, type)

Projects any nR x nC matrix in the equivariant subspace

Performs the integration

`` X1 = int{g in G} dg rhoR.image(g) * X * rhoC.inverseImage(g) ``

Note that there are little benefits usually to use the 'double/sparse' type compared to 'double'.

Raises

An error if type is 'exact' and isExact is false.

Parameters
  • X (double(*,*) or cyclotomic (*,*), may be sparse) – Matrix to project

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • X1 (double(*,*) or cyclotomic (*,*)) – Projected matrix

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned X1 to the invariant subspace in Frobenius norm

sample(self, type)

Returns an approximate sample from this equivariant space along with estimated numerical error

Raises

An error if type is 'exact' and isExact is false.

Parameters

type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • X (double(*,*)) – A sample from this equivariant space

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned X to the invariant subspace in Frobenius norm

subEquivariant(self, subR, subC, varargin)

Constructs a invariant subspace of an equivariant space

Parameters
  • subC (replab.SubRep) – A subrepresentation of self.repC

  • subR (replab.SubRep) – A subrepresentation of self.repR

  • special – Special structure if applicable, see Equivariant, default: ‘’

  • type – Whether to obtain an exact equivariant space, default ‘double’ (‘double’ and ‘double/sparse’ are equivalent)

Kwtype special

charstring, optional

Kwtype type

‘exact’, ‘double’ or ‘double/sparse’, optional

static validateArguments(repR, repC, special)

Validate and complete arguments used to define an equivariant space

When the special argument is specified, one can omit either repR or repC and it will be completed automatically.

Parameters
  • repR (Rep or []) – Representation acting on rows

  • repC (Rep or []) – Representation acting on columns

  • special (charstring or '') – One of the equivariant types, see Equivariant

Returns

  • repR (Rep) – Row representation

  • repC (Rep) – Column representation

SubEquivariant

class replab.SubEquivariant(parent, repR, repC, special)

Bases: replab.Equivariant

Describes a subspace of an equivariant space induced by subspaces of source/target representations

This equivariant space relates to a parent equivariant space as follows.

  • self.repR must be a SubRep of self.parent.repR,

  • self.repC must be a SubRep of self.parent.repC.

Own members

parent

Parent equivariant space

Type

Equivariant

static make(repR, repC, varargin)

Returns the space of equivariant linear maps between two subrepresentations

The equivariant vector space contains the matrices X such that

repC.image(g) * X = X * repR.image(g)

Parameters
  • repR (replab.SubRep) – Representation on the target/row space

  • repC (replab.SubRep) – Representation on the source/column space

  • special – Special structure if applicable, see Equivariant, default: ‘’

  • type – Whether to obtain an exact equivariant space, default ‘double’

  • parent – Equivariant space of repR.parent and repC.parent, default: []

Kwtype special

charstring, optional

Kwtype type

‘exact’, ‘double’ or ‘double/sparse’, optional

Kwtype parent

Equivariant, optional

Returns

The equivariant vector space

Return type

replab.SubEquivariant

projectFromParent(self, X, type)

Projects the given matrix, written in the parent representation space, into the equivariant space of the subrepresentation

Raises

An error if type is 'exact' and isExact is false.

Parameters
  • X (double(*,*) or cyclotomic (*,*), may be sparse) – Matrix to project

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • X1 (double(*,*) or cyclotomic (*,*)) – Projected matrix

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned X1 to the invariant subspace in Frobenius norm

IsotypicEquivariant

class replab.IsotypicEquivariant(parent, repR, repC, special, R_internal, divisionAlgebraName)

Bases: replab.SubEquivariant

Equivariant space between two harmonized isotypic components containing equivalent representations

We consider two cases.

If the two isotypic components contain equivalent irreducible representations, then isZero is false and the matrices in this equivariant space have the following form:

:math:` X = sum_i M(:,:,i) otimes R(:,:,i) otimes A(:,:,i) `

where

  • \(M\) represents the multiplicity space,

  • \(R\) is a constant matrix representing the representation space,

  • \(A\) is a constant matrix encoding the division algebra.

The division algebra matrix \(A\) is given by:

  • A = 1 if divisionAlgebraName is ''

  • A(:,:,1) = [1 0; 0 1], A(:,:,2) = [0 -1; 1 0] if divisionAlgebraName is 'C->R'

  • A(:,:,1) = [1 0 0 0; 0 1 0 0; 0 0 1 0; 0 0 0 1], A(:,:,2) = [0 -1 0 0; 1 0 0 0; 0 0 0 1; 0 0 -1 0], A(:,:,3) = [0 0 1 0; 0 0 0 1; -1 0 0 0; 0 -1 0 0], A(:,:,4) = [0 0 0 -1; 0 0 1 0; 0 -1 0 0; 1 0 0 0] if divisionAlgebraName is 'H->R:equivariant'.

If the two isotypic components contain inequivalent irreps, then isZero is true. divisionAlgebraName is '0', and R has size [m n 0].

Example

>>> S3 = replab.S(3);
>>> natDec = S3.naturalRep.decomposition;
>>> std1 = natDec.component(2);
>>> nat2Dec = S3.naturalRep.tensorPower(2).decomposition;
>>> std2 = nat2Dec.component(3);
>>> E = std1.isotypicEquivariantFrom(std2);
>>> X = E.sample;
>>> M = E.projectAndFactor(X);
>>> size(M, 3)
    1
>>> R = E.R('double');
>>> A = E.A('double');
>>> tol = 1e-10;
>>> norm(kron(M(:,:,1), kron(R(:,:,1), A(:,:,1))) - X, 'fro') <= tol
    1
>>> norm(E.reconstruct(M) - X, 'fro') <= tol
    1

Own members

R_internal

Representation space basis

Type

double(*,*,*) or cyclotomic (*,*,*)

divisionAlgebraName

Division algebra

Type

‘’, ‘C->R’, ‘H->R:equivariant’

A(self, type)

Returns the division algebra basis

Parameters

type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

The division algebra basis

Return type

double(*,*) or cyclotomic (*,*)

R(self, type)

Returns the representation space basis

Parameters

type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

The representation space basis

Return type

double(*,*,*) or cyclotomic (*,*,*)

divisionAlgebraBlockSize(self)

Returns the dimension of the matrix block encoding the division algebra

As a special case, if the block is zero, the returned value is one.

Returns

Dimension of the matrix block

Return type

integer

divisionAlgebraDimension(self)

Returns the dimension of the division algebra encoded in this block

As a special case, if the block is zero, the returned value is zero.

Returns

Dimension of the division algebra

Return type

integer

isZero(self)

Returns whether this equivariant space contains only the zero matrix

This happens when the isotypic components repR and repC correspond to inequivalent irreducible representations

Example

>>> S3 = replab.S(3);
>>> rep = S3.naturalRep;
>>> triv = rep.decomposition.component(1);
>>> std = rep.decomposition.component(2);
>>> E = triv.isotypicEquivariantFrom(std);
>>> E.isZero
    1
Returns

True if the equivariant space is trivial

Return type

logical

static kronFirstFactor(X, X2)

Factors a Kronecker product, second step

Assuming that X is approximately X = kron(X1, X2), with X2 known, this estimates the first factor X1 and returns it.

Parameters
  • X (double(*,*)) – Matrix to factor

  • X2 (double(*,*)) – Second factor

Returns

An estimation of the first factor in the Kronecker product

Return type

double(*,*)

static kronSecondFactor(X, nRows2, nCols2)

Factors a Kronecker product, first step

Assuming that X is approximately X = kron(X1, X2), with X2 a matrix of size nRows2 x nCols2, this estimates the factor X2 and returns it.

Parameters
  • X (double(*,*)) – Matrix to factor

  • nRows2 (integer) – Number of rows of the second factor

  • nCols2 (integer) – Number of columns of the second factor

Returns

An estimation of the second factor in the Kronecker product

Return type

double(*,*)

static make(repR, repC, varargin)

Returns the space of equivariant linear maps between two isotypic components

The equivariant vector space contains the matrices X such that

repC.image(g) * X = X * repR.image(g)

Parameters
  • repR (replab.Isotypic) – Isotypic component on the target/row space

  • repC (replab.Isotypic) – Isotypic component on the source/column space

  • special – Special structure, see help on Equivariant

  • type – Whether to obtain an exact equivariant space (‘double’ and ‘double/sparse’ are equivalent)

  • parent – Equivariant space from repC.parent to repR.parent

  • parentSample – A generic sample from the parent space

Kwtype special

charstring

Kwtype type

‘exact’, ‘double’ or ‘double/sparse’

Kwtype parent

Equivariant, optional

Kwtype parentSample

double(*,*), optional

Returns

The equivariant vector space

Return type

replab.IsotypicEquivariant

makeSdpvar(self)

Returns a parameterization of the multiplicity space using YALMIP variables

Returns

Parameterization

Return type

sdpvar(*,*,*)

projectAndFactor(self, X, type)

Projects the given matrix in the commutant algebra and factors it

It returns the decomposition of the projection.

Parameters
  • X (double(*,*) or cyclotomic (*,*), may be sparse) – Matrix in the isotypic component space to project

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • M (double(*,*,*) or cyclotomic (*,*,*)) – The part containing the degrees of freedom of the commutant algebra

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned projection to the invariant subspace in Frobenius norm

projectAndFactorFromParent(self, parentX, type)

Projects the given matrix in the parent representation space into the commutant algebra and factors it

It returns the decomposition of the projection.

Parameters
  • parentX (double(*,*) or cyclotomic (*,*), may be sparse) – Matrix in the parent representation space to project

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • M (double(*,*,*) or cyclotomic (*,*,*)) – The part containing the degrees of freedom of the commutant algebra

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned projection to the invariant subspace in Frobenius norm

reconstruct(self, M, type)

Reconstructs an equivariant space element returned by a factorization method

Parameters
  • M (double(*,*,*) or cyclotomic(*,*,*)) – Factorized element

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

Equivariant matrix

Return type

double(*,*) or cyclotomic(*,*)

IrreducibleEquivariant

class replab.IrreducibleEquivariant(parent, repR, repC, special, blocks)

Bases: replab.SubEquivariant

Equivariant space between two harmonized irreducible decompositions

Matrices in this equivariant space are composed of blocks with the following form:

\(X_{ij} = \sum_k M_{ijk} \otimes R_{ijk} \otimes A_{ijk\)} (or, most often, \(X_{ij}\) is identically zero)

where

  • \(M_{ijk}\) represents the multiplicity space,

  • \(R_{ijk}\) is a constant matrix representing the representation space,

  • \(A_{ijk}\) is a constant matrix encoding the division algebra.

Own members

blocks

Isotypic equivariant spaces

Type

cell(*,*) of IsotypicEquivariant

nonZeroBlock

True when the corresponding block can be nonzero

Type

logical(*,*)

colBlockSize(self)

Column sizes of blocks

Returns

Sizes

Return type

integer(1,*)

static make(repR, repC, varargin)

Returns the space of equivariant linear maps between two harmonized irreducible decompositions

The equivariant vector space contains the matrices X such that

repC.image(g) * X = X * repR.image(g)

Parameters
  • repR (replab.Irreducible) – Harmonized irreducible decomposition acting on the target/row space

  • repC (replab.Irreducible) – Harmonized irreducible decompositoin acting on the source/column space

  • special – Special structure, see help on Equivariant

  • type – Whether to obtain an exact equivariant space (‘double’ and ‘double/sparse’ are equivalent)

  • parent – Equivariant space from repC.parent to repR.parent

  • parentSample – A generic sample from the parent space

Kwtype special

charstring

Kwtype type

‘exact’, ‘double’ or ‘double/sparse’

Kwtype parent

Equivariant, optional

Kwtype parentSample

double(*,*), optional

Returns

The equivariant vector space

Return type

replab.IrreducibleEquivariant

makeSdpvarBlocks(self)

Returns a parameterization of the equivariant subspace using YALMIP variables

Returns

Parameterization using blocks

Return type

cell(1,*) of sdpvar(*,*,*)

projectAndFactor(self, X, type)

Projects the given matrix in the commutant algebra and factors it

It returns the decomposition of the projection.

Parameters
  • X (double(*,*) or cyclotomic (*,*), may be sparse) – Matrix in the isotypic component space to project

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • M (cell(*,*) of ([], or cell(1,*) of (double(*,*) or cyclotomic (*,*)))) – The part containing the degrees of freedom of the commutant algebra

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned projection to the invariant subspace in Frobenius norm

projectAndFactorFromParent(self, X, type)

Projects the given matrix in the parent representation space into the commutant algebra and factors it

It returns the decomposition of the projection.

Parameters
  • X (double(*,*) or cyclotomic (*,*), may be sparse) – Matrix in the parent representation space to project

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

  • M (cell(*,*) of ([], or cell(1,*) of (double(*,*) or cyclotomic (*,*)))) – The part containing the degrees of freedom of the commutant algebra

  • err (double) – Estimation of the numerical error, expressed as the distance of the returned projection to the invariant subspace in Frobenius norm

reconstruct(self, M, type)

Reconstructs an equivariant space element returned by a factorization method

Parameters
  • M (cell(*,*) of cell(1,*) of double(*,*) or cyclotomic(*,*)) – Factorized element

  • type ('double', 'double/sparse' or 'exact', optional) – Type of the returned value, default: ‘double’

Returns

Equivariant matrix

Return type

double(*,*) or cyclotomic(*,*)

rowBlockSize(self)

Row sizes of blocks

Returns

Sizes

Return type

integer(1,*)