Grassmann Operations

UQpy supports several operations on the Grassmann manifold. These operations are defined as methods of the GrassmannOperations class described herein.

class GrassmannOperations(grassmann_points, kernel=<UQpy.utilities.kernels.grassmannian_kernels.ProjectionKernel.ProjectionKernel object>, p='max', optimization_method='GradientDescent', distance=<UQpy.utilities.distances.grassmannian_distances.GeodesicDistance.GeodesicDistance object>)[source]

The GrassmannOperations class can used in two ways. In the first case, the user can invoke the initializer by providing all the required input data. The class will then automatically calculate the kernel matrix, distance matrix, karcher mean, and frechet variance of the input Grassmann points. Alternatively, the user may invoke each of the static methods individually and calculate only the required quantities, without having to instantiate a new object.

Parameters:
  • grassmann_points (Union[list[ndarray], list[GrassmannPoint]]) – Data points projected on the Grassmann manifold

  • kernel (GrassmannianKernel) – Kernel to be used to evaluate the kernel matrix of the given Grassmann points

  • p (Union[int, str]) – Rank of the Grassmann projected points.

  • optimization_method (str) –

    String that defines the optimization method for computation of the Karcher mean.

    Options: “GradientDescent”, “StochasticGradientDescent”

  • distance (GrassmannianDistance) – Distance measure to be used for the optimization in computation of the Karcher mean and Frechet variance.

Methods

The following sections introduce the methods available in the GrassmannOperations class with a brief introduction to their theory and their specifications in UQpy.

Exponential Map

For point \(\mathbf{X}\) on the Grassmann manifold \(\mathcal{G}(p,n)\), we can define the tangent space \(\mathcal{T}_{\mathbf{X}, \mathcal{G}(p,n)}\) as the set of all matrices \(\mathbf{\Gamma}\) such that

\[\mathcal{T}_{\mathbf{X}, \mathcal{G}(p,n)} = \{\mathbf{\Gamma} \in \mathbb{R}^{n \times p} :\mathbf{\Gamma}^T\mathbf{X}=\mathbf{0}\}\]

Consider two points \(\mathbf{X}\) and \(\mathbf{Y}\) on \(\mathcal{G}(p, n)\) and the tangent space \(\mathcal{T}_{\mathbf{X}, \mathcal{G}(p,n)}\) at \(\mathbf{X}\). The geodesic, which refers to the shortest curve on the manifold connecting the two points, is defined as

\[\mathbf{\Gamma} = \mathbf{U}\mathbf{S}\mathbf{V}^T\]
\[\Phi(t)=\mathrm{span}\left[\left(\mathbf{X}\mathbf{V}\mathrm{cos}(t\mathbf{S})+\mathbf{U}\mathrm{sin}(t\mathbf{S})\right)\mathbf{V}^T\right]\]

where \(\mathbf{\Phi}(0)=\mathbf{X}\) and \(\mathbf{\Phi}(1)=\mathbf{Y}\).

The exponential map, denoted as \(\mathrm{exp}_{\mathbf{X}}(\mathbf{\Gamma}): \mathcal{T}_{\mathbf{X}, \mathcal{G}(p,n)} \rightarrow \mathcal{G}(p, n)\), maps a tangent vector \(\mathbf{\Gamma} \in \mathcal{T}_{\mathbf{X}, \mathcal{G}(p,n)}\) to the endpoint \(\mathbf{Y}=\Phi(1)\) of the unique geodesic \(\Phi\) that emanates from \(\mathbf{X}\) in the direction \(\mathbf{\Gamma}\).

\[\mathrm{exp}_{\mathbf{X}}(\mathbf{\Gamma})\equiv\mathbf{Y}\]
\[\mathbf{Y} = \mathrm{exp}_{\mathbf{X}}(\mathbf{U}\mathbf{S}\mathbf{V}^T) = \mathbf{X}\mathbf{V}\mathrm{cos}\left(\mathbf{S}\right)\mathbf{Q}^T+\mathbf{U}\mathrm{sin}\left(\mathbf{S}\right)\mathbf{Q}^T\]

The exponential map is implemented in UQpy through the static exp_map() method.

To use the exp_map() method, one needs to import the GrassmannOperations class from the UQpy.dimension_reduction.grassmann_manifold module as follows:

>>> from UQpy.dimension_reduction.grassmann_manifold import GrassmannOperations

Since exp_map() is a static method, it does not require instantiation of the GrassmannOperations class.

static GrassmannOperations.exp_map(tangent_points, reference_point)[source]

Maps tangent vector(s) to the endpoint of a unique geodesic.

Parameters:
Return type:

list[GrassmannPoint]

Returns:

Point(s) on the Grassmann manifold.

Logarithmic Map

The logarithmic map, denoted as \(\mathrm{log}_\mathbf{X}(\mathbf{Y}):\mathcal{G}(p, n) \rightarrow \mathcal{T}_{\mathbf{X}\mathcal{G}(p,n)}\) maps the endpoint \(\mathbf{Y}=\Phi(1)\) of the unique geodesic \(\Phi\) that emanates from \(\mathbf{X}\) to a tangent vector \(\mathbf{\Gamma} \in \mathcal{T}_{\mathbf{X}, \mathcal{G}(p,n)}\).

\[\mathrm{log}_\mathbf{X}(\mathbf{Y})\equiv \mathbf{\Gamma} = \mathbf{U}\mathrm{tan}^{-1}\left(\mathbf{S}\right)\mathbf{V}^T\]

where \(\mathbf{U}, \mathbf{S}, \mathbf{V}\) are define as in the exponential map above.

To use the log_map() method, one needs to import the GrassmannOperations class from the UQpy.dimension_reduction.grassmann_manifold module as follows:

>>> from UQpy.dimension_reduction.grassmann_manifold import GrassmannOperations

Since log_map() is a static method, it does not require instantiation of the GrassmannOperations class.

static GrassmannOperations.log_map(grassmann_points, reference_point)[source]

Maps the endpoint of the unique geodesic from \(\mathbf{X}\) to tangent vector(s)

Parameters:
Return type:

list[ndarray]

Returns:

Point(s) on the tangent space.

Karcher mean

For a set of points \(\{\mathbf{X}_i\}_{i=1}^N\) on \(\mathcal{G}(p,n)\) the Karcher mean is defined as the solution of the following minimization problem:

\[\mu = \arg_{\mathbf{Y}} \mathrm{min}\left(\frac{1}{N}\sum_{i=1}^N d(\mathbf{X}_i, \mathbf{Y})^2\right)\]

where \(d(\cdot)\) is a Grassmann distance measure and \(\mathbf{Y}\) is a reference point on \(\mathcal{G}(p,n)\).

To use the karcher_mean() method, one needs to import the GrassmannOperations class from the UQpy.dimension_reduction.grassmann_manifold module as follows:

>>> from UQpy.dimension_reduction.grassmann_manifold import GrassmannOperations

Since karcher_mean() is a static method, it does not require instantiation of the GrassmannOperations class.

static GrassmannOperations.karcher_mean(grassmann_points, optimization_method, distance, acceleration=False, tolerance=0.001, maximum_iterations=1000)[source]
Parameters:
  • maximum_iterations (int) – Maximum number of iterations performed by the optimization algorithm.

  • tolerance (float) – Tolerance used as the convergence criterion of the optimization.

  • acceleration (bool) – Boolean flag used in combination with GradientDescent optimization method which activates the Nesterov acceleration scheme

  • grassmann_points (Union[list[ndarray], list[GrassmannPoint]]) – Point(s) on the Grassmann manifold.

  • optimization_method (str) –

    String that defines the optimization method.

    Options: “GradientDescent”, “StochasticGradientDescent”

  • distance (GrassmannianDistance) – Distance measure to be used for the optimization.

Return type:

GrassmannPoint

UQpy offers two methods for solving this optimization, the GradientDescent and the StochasticGradientDescent. Both are implemented as private methods of the GrassmannOperations class.

Frechet variance

For a set of points \(\{\mathbf{X}_i\}_{i=1}^N\) on \(\mathcal{G}(p,n)\) the Frechet variance is defined as:

\[\sigma_{f}^2 = \frac{1}{N}\sum_{i=1}^N d(\mathbf{X}_i, \mu)^2\]

where \(d(\cdot)\) is a Grassmann distance measure and \(\mu\) is the Karcher mean of set of points \(\{\mathbf{X}_i\}_{i=1}^N\) on \(\mathcal{G}(p,n)\).

To use the frechet_variance() method, one needs to import the GrassmannOperations class from the UQpy.dimension_reduction.grassmann_manifold module as follows:

>>> from UQpy.dimension_reduction.grassmann_manifold import GrassmannOperations

Since frechet_variance() is a static method, it does not require instantiation of the GrassmannOperations class.

static GrassmannOperations.frechet_variance(grassmann_points, reference_point, distance)[source]

Variance measure of the Grassman distance in relation to the Karcher mean.

Parameters:
Return type:

float