geomstats package#

Subpackages#

Submodules#

geomstats.algebra_utils module#

Utility module of reusable algebra routines.

geomstats.algebra_utils.flip_determinant(matrix, det)[source]#

Change sign of the determinant if it is negative.

For a batch of matrices, multiply the matrices which have negative determinant by a diagonal matrix :math:`diag(1,…,1,-1) from the right. This changes the sign of the last column of the matrix.

Parameters:
  • matrix (array-like, shape=[…,n ,m]) – Matrix to transform.

  • det (array-like, shape=[…]) – Determinant of matrix, or any other scalar to use as threshold to determine whether to change the sign of the last column of matrix.

Returns:

matrix_flipped (array-like, shape=[…, n, m]) – Matrix with the sign of last column changed if det < 0.

geomstats.algebra_utils.from_vector_to_diagonal_matrix(vector, num_diag=0)[source]#

Create diagonal matrices from rows of a matrix.

Parameters:
  • vector (array-like, shape=[m, n])

  • num_diag (int) – number of diagonal in result matrix. If 0, the result matrix is a diagonal matrix; if positive, the result matrix has an upper-right non-zero diagonal; if negative, the result matrix has a lower-left non-zero diagonal. Optional, Default: 0.

Returns:

diagonals (array-like, shape=[m, n, n]) – 3-dimensional array where the i-th n-by-n array diagonals[i, :, :] is a diagonal matrix containing the i-th row of vector.

geomstats.algebra_utils.rotate_points(points, end_point)[source]#

Apply to points the rotation from north_pole to end_point.

A QR decomposition is used to find the rotation that maps the north pole (1, 0,…,0) to the end_point, then this rotation is applied to the input points.

Parameters:
  • points (array-like, shape=[…, n]) – Points to rotate.

  • end_point (array-like, shape=[n, ]) – Point to parametrise the rotation.

Returns:

rotated_points (array-like, shape=[…, n]) – Points after the rotation.

geomstats.algebra_utils.taylor_exp_even_func(point, taylor_function, order=5, tol=1e-06)[source]#

Taylor Approximation of an even function around zero.

Parameters:
  • point (array-like) – Argument of the function to approximate.

  • taylor_function (dict with following keys) –

    functioncallable

    Even function to approximate around zero.

    coefficientslist

    Taylor coefficients of even order at zero.

  • order (int) – Order of the Taylor approximation. Optional, Default: 5.

  • tol (float) – Threshold to use the approximation instead of the function’s value. Where abs(point) <= tol, the approximation is returned.

Returns:

function_value (array-like) – Value of the function at point.

geomstats.errors module#

Checks and associated errors.

exception geomstats.errors.ShapeError[source]#

Bases: ValueError

Raised when there is an incompatibility between shapes.

geomstats.errors.check_belongs(point, manifold, atol=1e-12)[source]#

Raise an error if point does not belong to the input manifold.

Parameters:
  • point (array-like) – Point to be tested.

  • manifold (Manifold) – Manifold to which the point should belong.

  • manifold_name (string) – Name of the manifold for the error message.

geomstats.errors.check_integer(n, n_name)[source]#

Raise an error if n is not a > 0 integer.

Parameters:
  • n (unspecified) – Parameter to be tested.

  • n_name (string) – Name of the parameter.

geomstats.errors.check_parameter_accepted_values(param, param_name, accepted_values)[source]#

Raise an error if parameter does not belong to a set of values.

Parameters:
  • param (unspecified) – Parameter to be tested.

  • param_name (string) – Name of the parameter.

  • accepted_values (list) – Accepted values that the parameter can take.

geomstats.errors.check_point_shape(point, manifold, suppress_error=False)[source]#

Check if the shape of point does not match the shape of a manifold or metric.

If the final elements of the shape of point do not match the shape of manifold (which may be any object with a shape attribute, such as a Riemannian metric) then point cannot be an array of points on the manifold (or similar) and a ValueError is raised. The error can be suppressed by setting suppress_error to True.

Parameters:
  • point (array-like) – The point to check the shape of.

  • manifold ({Manifold, RiemannianMetric}) – The object to check the point against

  • suppress_error (bool) – Whether to suppress the ShapeError if the shapes do not match. Optional, default is False.

Returns:

shapes_match (bool) – Whether the shape of the point matches the shape of the manifold or metric.

Raises:

ValueError – If the final dimensions of point are not equal to the final dimensions of manifold.

geomstats.errors.check_positive(param, param_name)[source]#

Raise an error if param is not a > 0 number.

Parameters:
  • param (unspecified) – Parameter to be tested.

  • param_name (string) – Name of the parameter.

geomstats.exceptions module#

Geomstats custom exceptions.

exception geomstats.exceptions.AutodiffNotImplementedError[source]#

Bases: RuntimeError

Raised when autodiff is not implemented.

exception geomstats.exceptions.NotPartialOrder[source]#

Bases: Exception

Raise an exception when less equal is not true.

geomstats.integrator module#

Integrator functions used when no closed forms are available.

Lead author: Nicolas Guigui.

These are designed for first order systems of ODEs written as a spatial variable x and a time variable t:

\[\frac{dx}{dt} = force(x, t)\]

where \(x\) is called the state variable. It may represent many variables by stacking arrays, e.g. position and velocity in a geodesic equation.

geomstats.integrator.euler_step(force, state, time, dt)[source]#

Compute one step of the euler approximation.

Parameters:
  • force (callable) – Vector field that is being integrated.

  • state (array-like, shape=[…, n, dim]) – State at time t, corresponds to position and velocity variables at time t.

  • time (float) – Time variable.

  • dt (float) – Time-step in the integration.

Returns:

point_new (array-like, shape=[…, n, dim]) – Variables at time t + dt.

geomstats.integrator.integrate(function, initial_state, end_time=1.0, n_steps=10, step='euler')[source]#

Compute the flow under the vector field using symplectic euler.

Integration function to compute flows of vector fields on a regular grid between 0 and a finite time from an initial state.

Parameters:
  • function (callable) – Vector field to integrate.

  • initial_state (tuple of arrays) – Initial position and speed.

  • end_time (float) – Final integration time. Optional, default : 1.

  • n_steps (int) – Number of integration steps to use. Optional, default : 10.

  • step (str, {‘euler’, ‘rk4’, ‘group_rk2’, ‘group_rk4’}) – Numerical scheme to use for elementary integration steps. Optional, default : ‘euler’.

Returns:

final_state (tuple) – sequences of solutions every end_time / n_steps. The shape of each element of the sequence is the same as the vectors passed in initial_state.

geomstats.integrator.leapfrog_step(force, state, time, dt)[source]#

Compute one step of the leapfrog approximation.

Parameters:
  • state (array-like, shape=[…, 2, dim]) – State at time t, corresponds to position and velocity variables at time t.

  • force (callable) – Vector field that is being integrated.

  • time (float) – Time variable.

  • dt (float) – Time-step in the integration.

Returns:

state_new (array-like, shape=[…, 2, dim]) – State at time t + dt, corresponds to position and velocity variables at time t + dt.

See also

https

//en.wikipedia.org/wiki/Leapfrog_integration

geomstats.integrator.rk2_step(force, state, time, dt)[source]#

Compute one step of the rk2 approximation.

Parameters:
  • force (callable) – Vector field that is being integrated.

  • state (array-like, shape=[…, n, dim]) – State at time t, corresponds to position and velocity variables at time t.

  • time (float) – Time variable.

  • dt (float) – Time-step in the integration.

Returns:

point_new (array-like, shape=[…, n, dim]) – Variables at time t + dt.

See also

https

//en.wikipedia.org/wiki/Runge–Kutta_methods

geomstats.integrator.rk4_step(force, state, time, dt)[source]#

Compute one step of the rk4 approximation.

Parameters:
  • force (callable) – Vector field that is being integrated.

  • state (array-like, shape=[…, n, dim]) – State at time t, corresponds to position and velocity variables at time t.

  • time (float) – Time variable.

  • dt (float) – Time-step in the integration.

Returns:

point_new (array-like, shape=[…, n, dim]) – Variables at time t + dt.

See also

https

//en.wikipedia.org/wiki/Runge–Kutta_methods

geomstats.integrator.symplectic_euler_step(force, state, time, dt)[source]#

Compute one step of the symplectic euler approximation.

Parameters:
  • state (array-like, shape=[…, 2, dim]) – State at time t, corresponds to position and velocity variables at time t.

  • force (callable) – Vector field that is being integrated.

  • time (float) – Time variable.

  • dt (float) – Time-step in the integration.

Returns:

  • point_new (array-like, shape=[…, 1, dim]) – Position variable at time t + dt.

  • vector_new (array-like, shape=[…, 1, dim]) – Velocity variable at time t + dt.

geomstats.vectorization module#

Decorator to handle vectorization.

This abstracts the backend type.

geomstats.vectorization.broadcast_to_multibatch(batch_shape_a, batch_shape_b, array_a, *array_b)[source]#

Broadcast to multibatch.

Gives to both arrays batch shape batch_shape_b + batch_shape_a.

Does nothing if one of the batch shapes is empty.

Parameters:
  • batch_shape_a (tuple) – Batch shape of array_a.

  • batch_shape_b (tuple) – Batch shape of array_b.

  • array_a (array)

  • array_b (array)

geomstats.vectorization.check_is_batch(point_ndim, *point)[source]#

Check if inputs are batch.

Parameters:
  • point_ndim (int) – Point number of array dimensions.

  • point (array-like) – Point belonging to the space.

Returns:

is_batch (bool) – Returns True if point contains several points.

geomstats.vectorization.get_batch_shape(point_ndim, *point)[source]#

Get batch shape.

Parameters:
  • point_ndim (int) – Point number of array dimensions.

  • point (array-like or None) – Point belonging to the space.

Returns:

batch_shape (tuple) – Returns the shape related with batch. () if only one point.

geomstats.vectorization.get_n_points(point_ndim, *point)[source]#

Compute the number of points.

Parameters:
  • point_ndim (int) – Point number of array dimensions.

  • point (array-like) – Point belonging to the space.

Returns:

n_points (int) – Number of points.

geomstats.vectorization.repeat_out(point_ndim, out, *point, out_shape=())[source]#

Repeat out shape after finding batch shape.

Parameters:
  • point_ndim (int) – Point number of array dimensions.

  • out (array-like) – Output to be repeated

  • point (array-like or None) – Point belonging to the space.

  • out_shape (tuple) – Indicates out shape for no batch computations.

Returns:

out (array-like) – If no batch, then input is returned. Otherwise it is broadcasted.

geomstats.vectorization.repeat_out_multiple_ndim(out, point_ndim_1, points_1, point_ndim_2, points_2, out_ndim=0)[source]#

Repeat out after finding batch shape.

Differs from repeat_out by accepting two sets of point_ndim arrays.

Parameters:
  • out (array-like) – Output to be repeated

  • point_ndim_1 (int) – Point number of array dimensions.

  • points_1 (tuple[array-like or None]) – Arrays of dimension point_ndim_1 or higher.

  • point_ndim_2 (int) – Point number of array dimensions.

  • points_2 (tuple[array-like or None]) – Arrays of dimension point_ndim_2 or higher.

  • out_ndim (int) – Out number of array dimensions.

Returns:

out (array-like) – If no batch, then input is returned. Otherwise it is broadcasted.

geomstats.vectorization.repeat_point(point, n_reps=2, expand=False)[source]#

Repeat point.

Parameters:
  • point (array-like) – Point of a space.

  • n_reps (int) – Number of times the point should be repeated.

  • expand (bool) – Repeat even if n_reps == 1.

Returns:

rep_point (array-like) – point repeated n_reps times.

Module contents#

Import main modules.