pyphysim.subspace package

Submodules

pyphysim.subspace.metrics module

Implement several metrics for subspaces.

pyphysim.subspace.metrics.calc_chordal_distance(matrix1: numpy.ndarray, matrix2: numpy.ndarray)float[source]

Calculates the chordal distance between the two matrices

Parameters
  • matrix1 (np.ndarray) – A 2D numpy array.

  • matrix2 (np.ndarray) – A 2D numpy array.

Returns

chord_dist – The chordal distance.

Return type

float

Notes

Same as calc_chordal_distance_2(), but implemented differently.

Examples

>>> A = np.arange(1, 9.)
>>> A.shape = (4, 2)
>>> B = np.array([[1.2, 2.1], [2.9, 4.3], [5.2, 6.1], [6.8, 8.1]])
>>> print(round(calc_chordal_distance(A, B), 8))
0.47386786
pyphysim.subspace.metrics.calc_chordal_distance_2(matrix1: numpy.ndarray, matrix2: numpy.ndarray)float[source]

Calculates the chordal distance between the two matrices

Parameters
  • matrix1 (np.ndarray) – A 2D numpy array.

  • matrix2 (np.ndarray) – A 2D numpy array.

Returns

chord_dist – The chordal distance.

Return type

float

Notes

Same as calc_chordal_distance(), but implemented differently.

Examples

>>> A = np.arange(1, 9.)
>>> A.shape = (4, 2)
>>> B = np.array([[1.2, 2.1], [2.9, 4.3], [5.2, 6.1], [6.8, 8.1]])
>>> print(round(calc_chordal_distance_2(A, B), 8))
0.47386786
pyphysim.subspace.metrics.calc_chordal_distance_from_principal_angles(principalAngles: numpy.ndarray)float[source]

Calculates the chordal distance from the principal angles.

It is given by the square root of the sum of the squares of the sin of the principal angles.

Parameters

principalAngles (np.ndarray) – Numpy array with the principal angles. This is a 1D numpy array.

Returns

chord_dist – The chordal distance.

Return type

float

Examples

>>> A = np.arange(1, 9.)
>>> A.shape = (4, 2)
>>> B = np.array([[1.2, 2.1], [2.9, 4.3], [5.2, 6.1], [6.8, 8.1]])
>>> princ_angles = calc_principal_angles(A, B)
>>> print(round(calc_chordal_distance_from_principal_angles(princ_angles), 8))
0.47386786
pyphysim.subspace.metrics.calc_principal_angles(matrix1: numpy.ndarray, matrix2: numpy.ndarray)numpy.ndarray[source]

Calculates the principal angles between matrix1 and matrix2.

Parameters
  • matrix1 (np.ndarray) – A 2D numpy array.

  • matrix2 (np.ndarray) – A 2D numpy array.

Returns

The principal angles between matrix1 and matrix2. This is a 1D numpy array.

Return type

np.ndarray

Examples

>>> A = np.array([[1, 2], [3, 4], [5, 6]])
>>> B = np.array([[1, 5], [3, 7], [5, -1]])
>>> print(calc_principal_angles(A, B))
[0.         0.54312217]

pyphysim.subspace.projections module

Module related to subspace projection.

class pyphysim.subspace.projections.Projection(A: numpy.ndarray)[source]

Bases: object

Class to calculate the projection, orthogonal projection and reflection of a given matrix in a Subspace S spanned by the columns of a matrix A.

The matrix A is provided in the constructor and after that the functions project, oProject and reflect can be called with M as an argument.

Parameters

A (np.ndarray) – The matrix whose columns form a basis for the projected subspace.

Examples

>>> A = np.array([[1 + 1j, 2 - 2j], [3 - 2j, 0], [-1 - 1j, 2 - 3j]])
>>> v = np.array([1, 2, 3])
>>> P = Projection(A)
>>> P.project(v)
array([1.69577465+0.87887324j, 1.33802817+0.41408451j,
       2.32957746-0.56901408j])
>>> P.oProject(v)
array([-0.69577465-0.87887324j,  0.66197183-0.41408451j,
        0.67042254+0.56901408j])
>>> P.reflect(v)
array([-2.3915493 -1.75774648j, -0.67605634-0.82816901j,
       -1.65915493+1.13802817j])
static calcOrthogonalProjectionMatrix(A: numpy.ndarray)numpy.ndarray[source]

Calculates the projection matrix that projects a vector (or a matrix) into the signal space orthogonal to the signal space spanned by the columns of M.

Parameters

A (np.ndarray) – A matrix whose columns form a basis for the “desired subspace”.

Returns

The projection matrix that can be used to project a vector or a matrix into the subspace orthogonal to the subspace spanned by the columns of A

Return type

np.ndarray

Examples

>>> A = np.array([[1, 2], [2, 2], [4, 3]])
>>> # Matrix that projects into the subspace orthogonal to the
>>> # subspace spanned by the columns of A
>>> oQ = calcOrthogonalProjectionMatrix(A)
>>> print(oQ)
[[ 0.12121212 -0.3030303   0.12121212]
 [-0.3030303   0.75757576 -0.3030303 ]
 [ 0.12121212 -0.3030303   0.12121212]]
static calcProjectionMatrix(A: numpy.ndarray)numpy.ndarray[source]

Calculates the projection matrix that projects a vector (or a matrix) into the signal space spanned by the columns of A.

Parameters

A (np.ndarray) – A matrix whose columns form a basis for the desired subspace.

Returns

The projection matrix that can be used to project a vector or a matrix into the subspace spanned by the columns of A

Return type

np.ndarray

Examples

>>> A = np.array([[1 + 1j, 2 - 2j], [3 - 2j, 0],                           [-1 - 1j, 2 - 3j]])
>>> # Matrix that projects into the subspace spanned by the columns
>>> # of A
>>> Q = calcProjectionMatrix(A)
>>> np.allclose(Q.round(4), np.array(           [[ 0.5239+0.j, 0.0366+0.3296j, 0.3662+0.0732j],            [ 0.0366-0.3296j, 0.7690+0.j, -0.0789+0.2479j],            [ 0.3662-0.0732j, -0.0789-0.2479j, 0.7070-0.j]]))
True
oProject(M: numpy.ndarray)numpy.ndarray[source]

Project the matrix (or vector) M the subspace ORTHOGONAL to the subspace projected with project.

Parameters

M (np.ndarray) – The matrix to be projected.

Returns

The projection of M into the orthogonal subspace.

Return type

np.ndarray

project(M: numpy.ndarray)numpy.ndarray[source]

Project the matrix (or vector) M in the desired subspace.

Parameters

M (np.ndarray) – The matrix to be projected.

Returns

The projection of M into the desired subspace.

Return type

np.ndarray

reflect(M: numpy.ndarray)numpy.ndarray[source]

Find the reflection of the matrix in the subspace spanned by the columns of A

Parameters

M (np.ndarray) – The matrix to be projected.

Returns

The reflection of M in the subspace.

Return type

np.ndarray

pyphysim.subspace.projections.calcOrthogonalProjectionMatrix(A: numpy.ndarray)numpy.ndarray

Calculates the projection matrix that projects a vector (or a matrix) into the signal space orthogonal to the signal space spanned by the columns of M.

Parameters

A (np.ndarray) – A matrix whose columns form a basis for the “desired subspace”.

Returns

The projection matrix that can be used to project a vector or a matrix into the subspace orthogonal to the subspace spanned by the columns of A

Return type

np.ndarray

Examples

>>> A = np.array([[1, 2], [2, 2], [4, 3]])
>>> # Matrix that projects into the subspace orthogonal to the
>>> # subspace spanned by the columns of A
>>> oQ = calcOrthogonalProjectionMatrix(A)
>>> print(oQ)
[[ 0.12121212 -0.3030303   0.12121212]
 [-0.3030303   0.75757576 -0.3030303 ]
 [ 0.12121212 -0.3030303   0.12121212]]
pyphysim.subspace.projections.calcProjectionMatrix(A: numpy.ndarray)numpy.ndarray

Calculates the projection matrix that projects a vector (or a matrix) into the signal space spanned by the columns of A.

Parameters

A (np.ndarray) – A matrix whose columns form a basis for the desired subspace.

Returns

The projection matrix that can be used to project a vector or a matrix into the subspace spanned by the columns of A

Return type

np.ndarray

Examples

>>> A = np.array([[1 + 1j, 2 - 2j], [3 - 2j, 0],                           [-1 - 1j, 2 - 3j]])
>>> # Matrix that projects into the subspace spanned by the columns
>>> # of A
>>> Q = calcProjectionMatrix(A)
>>> np.allclose(Q.round(4), np.array(           [[ 0.5239+0.j, 0.0366+0.3296j, 0.3662+0.0732j],            [ 0.0366-0.3296j, 0.7690+0.j, -0.0789+0.2479j],            [ 0.3662-0.0732j, -0.0789-0.2479j, 0.7070-0.j]]))
True

Module contents

Subspace related stuff.

This includes projections and metrics for subspaces.