ddg.math.euclidean module

Utility functions for Euclidean space.

It includes parametrizations for circles/spheres in three-dimensional space and an implementation of Catmull–Rom splines.

ddg.math.euclidean.torus_parametrization(phi, theta, R, r)[source]

A simple parametrization of the torus.

Parameters:
phifloat

The value of the first parameter.

thetafloat

The value of the second parameter.

Rfloat

Radius of the larger circle.

rfloat

Radius of the smaller circle.

Returns:
np.ndarray of shape (3, 1)
ddg.math.euclidean.mobius_strip_parametrization(u, v, R)[source]

A simple parametrization of the moebius strip.

Parameters:
ufloat

The value of the first parameter.

vfloat

The value of the second parameter.

Rfloat

Radius of the midcircle of the moebius strip.

Returns:
np.ndarray of shape (3, 1)
ddg.math.euclidean.intersect_diags(v1, v2, v3, v4)[source]

For a quadrilateral(v1,v2,v3,v4) that lies in a plane, get the intersection point of the diagonals.

The four points are given as an argument in positive cyclic order.

v4-----v3
 |\   /|
 |  x  |
 |/   \|
v1-----v2
Parameters:
v1, v2, v3, v4numpy.ndarray of shape (n,)

Four points of quadrilateral vertices in positive cyclic order.

Returns:
numpy.ndarray of shape (n,)

Notes

Works for arbitrarily many dimensions as long as it is contained in a plane.

ddg.math.euclidean.intersect_edges(v1, v2, v3, v4)[source]

For a quadrilateral (v1,v2,v3,v4) that lies in a plane, get the intersection of the two opposite edges (v1,v2) and (v3,v4).

v4-----v3--
 |     |    \
 |     |     x
 |     |    /
v1-----v2--
Parameters:
v1, v2, v3, v4numpy.ndarray of shape (n,)

Four points of quadrilateral vertices in positive cyclic order.

Returns:
numpy.ndarray of shape (n,)

Notes

Works for arbitrarily many dimensions as long as it is contained in a plane.

ddg.math.euclidean.christoffel_dual_quad(v1, v2, v3, v4)[source]

Compute the Christoffel dual quadrilateral of a given quadrilateral.

Corresponding edges and non-corresponding diagonals are parallel. This determines the dual quadrilateral up to scaling.

Parameters:
v1, v2, v3, v4numpy.ndarray of shape (n,)

Coordinates of the given quadrilateral.

Returns:
V1, V2, V3, V4numpy.ndarray of shape (n,)

Coordinates of the dual quadrilateral.

ddg.math.euclidean.focal_points(v1, v2, v3, v4)[source]

For a quadrilateral (v1, v2, v3, v4), get the focal points of the quadrilateral.

These are the intersecting points of all edges.

   p2
   / \
  /   \
v4-----v3--
 |     |    \
 |     |     p1
 |     |    /
v1-----v2--
Parameters:
v1, v2, v3, v4numpy.ndarray of shape (n,)

Coordinates of the given quadrilateral.

Returns:
p1, p2numpy.ndarray of shape (n,)

Focal points of the given quadrilateral.

ddg.math.euclidean.diagonal_triangle(v1, v2, v3, v4)[source]

For a quadrilateral (v1, v2, v3, v4), get the diagonal triangle of the quadrilateral.

These points determine the triangle formed by the intersection of all diagonals.

   p3
   / \
  /   \
v4-----v3--
 |\   /|    \
 |  p1 |     p2
 |/   \|    /
v1-----v2--
Parameters:
v1, v2, v3, v4numpy.ndarray of shape (n,)

Coordinates of the given quadrilateral.

Returns:
p1, p2, p3numpy.ndarray of shape (n,)

Diagonal triangle of the given quadrilateral.

ddg.math.euclidean.circle_fct(t, center, radius=1, normals=None)[source]

Parametrization map of a circle with a given center, radius and list of normals. The dimensions of the normals must match the dimension of the center and both must be >= 2.

Parameters:
t: float

Parameter for the parametrization.

center: numpy.ndarray of shape (d,)

Center of the circle.

radius: float (default=1)

Radius of the circle.

normals: numpy.ndarray of shape (n,d) (default=None)

Normals of the circle that determine the orthogonal complement of the subspace containing the circle. If None are given the orthogonal complement will be determined by the dimension of the center and its normals consist of unit normal vectors. These vectors must be linearly independent(!)

Returns:
float

Value of a point on the circle at parameter t.

Raises:
ValueError

If the shapes of center and normals do not match. If the number of normals is not smaller than the dimension of the circle center minus 2.

ddg.math.euclidean.sphere_fct(u, v, center=None, radius=1)[source]

Parameterization map of a 2d sphere in 3d euclidean space with a given center and radius

Parameters:
u: float

First parameter for the parametrization.

v: float

Second parameter for the parametrization.

center: array_like of shape (3,), (default=None)

Center of the sphere.

radius: float (default=1)

Radius of the sphere.

Returns:
float

Value of the sphere parametrization for given parameters u and v.

Raises:
ValueError

If center is not 3-dimensional.

ddg.math.euclidean.circumcenter(v0, v1, v2)[source]

Computes the center of the circle through the three given points in two-dimensional or three-dimensional space.

Parameters:
v0, v1, v2array_like of shape (2,) or (3,)

Points in R^2 or R^3

Returns:
numpy.ndarray of shape (2,) or (3,)

Center of the circle through the points v0, v1, v2

ddg.math.euclidean.circle_through_three_points(v0, v1, v2, atol=None, rtol=None)[source]

Computes the center, radius and normal of the circle through the three given points in three-dimensional space.

Parameters:
v0, v1, v2array_like of shape (2,) or (3,)

Points in R^2 or R^3.

atol, rtolfloat (default=None)

If None, the global defaults will be used.

Returns:
numpy.ndarray of shape (2,) (3,), float, numpy.ndarray of shape (3,) or float

Center, radius, normal of the circle through the three given points, respectively. If points are 2D then center will have shape (2,) and normal will be float.

Raises:
ValueError

If the three points given lie on a common line

Notes

Since the cross product of difference vectors of the points is computed, their dimension must be three.

If 2d points are given instead, the third component is assumed to be zero and the normal is returned as a value (type float), such that the normal vector is given by [0, 0, normal]

ddg.math.euclidean.sphere_through_four_points(p0, p1, p2, p3)[source]

Computes the center and radius of sphere through the four given points in three-dimensional space.

Parameters:
p0, p1, p2, p3array_like of shape (3,)

Points in R^3.

Returns:
numpy.ndarray of shape (3,), float

Center, radius of the sphere respectively.

ddg.math.euclidean.rotation_2d(alpha)[source]

Returns 2x2 rotation matrix

      cos(alpha) | -sin(alpha)
A =   ------------------------
      sin(alpha) | cos(alpha)
Parameters:
alphafloat

Rotation angle.

Returns:
numpy.ndarray of shape (2, 2)
ddg.math.euclidean.rotation_from_to(source, target, homogeneous=True, atol=None, rtol=None)[source]

Returns a 3x3 rotation matrix taking the source vector to the target vector, or a 4x4 matrix of the type

      A | 0
T =   -----
      0 | 1
Parameters:
sourcearray_like of shape (3,)

Source vector.

targetarray_like of shape (3,)

Target vector.

homogeneous: bool, (default=True)

Determines whether 4x4 matrix or 3x3 matrix is returned.

atol, rtolfloat (default=None)

If None, the global defaults will be used.

Returns:
numpy.ndarray of shape either (3, 3) or (4, 4)

Rotation matrix taking source vector to target vector.

Raises:
ValueError

If at least one of the vectors is the zero vector (0, 0, 0).

Notes

If source and target are parallel the identety matirx is returned. If source and target are anti-parallel a correspoding 180 degree rotation is returned.

ddg.math.euclidean.scale_rotation_from_to(source, target, homogeneous=True, atol=None, rtol=None)[source]

Returns a 3x3 scale rotation matrix taking the source vector to the target vector, or a 4x4 matrix of the type

      A | 0
T =   -----
      0 | 1
Parameters:
sourcearray_like of shape (3,)

Source vector.

targetarray_like of shape (3,)

Target vector.

homogeneous: bool, (default=True)

Determines whether 4x4 matrix or 3x3 matrix is returned.

atol, rtolfloat (default=None)

If None, the global defaults will be used.

Returns:
numpy.ndarray of shape either (3, 3) or (4, 4)

Scale rotation matrix taking source vector to target vector.

Raises:
ValueError

If at least one of the vectors is the zero vector (0, 0, 0).

ddg.math.euclidean.rotation_angle_axis(rotation_vector, homogeneous=True)[source]

This function parameterizes rotation using axis–angle representation.

It returns either 4x4 matrix of the type

      A | 0
M =   -----
      0 | 1

where A represents the rotation matrix or A as a 3x3 rotation matrix.

Parameters:
rotation_vectorarray_like of shape (3,)

A 3-dimensional vector which is co-directional to the axis of rotation and whose norm gives the angle of rotation.

homogeneous: bool, (default=True)

Determines whether 4x4 matrix or 3x3 matrix is returned.

Returns:
numpy.ndarray of shape (3, 3) or (4, 4)
ddg.math.euclidean.scaleXYZ(xyz=(1, 1, 1), homogeneous=True)[source]

Returns either 4x4 matrix of type

      A | 0
M =   -----
      0 | 1

where A denotes diagonal scaling matrix, or A as 3x3 diagonal scaling matrix.

Parameters:
xyzarray_like of shape (3,), (default=(1, 1, 1))

Scaling factors in x, y and z directions.

homogeneous: bool, (default=True)

Determines whether 4x4 matrix or 3x3 matrix is returned.

Returns:
numpy.ndarray of shape (3, 3) or (4, 4)

Scaling matrix in either affine or homogeneous coordinates.

ddg.math.euclidean.translation_to(b=(0, 0, 0))[source]

Returns 4x4 matrix of type

      I | b
M =   -----
      0 | 1

where I denotes 3x3 identity matrix and b denotes 3x1 translation vector.

Parameters:
barray_like of shape (3,), (default=(0, 0, 0))

Translation vector.

Returns:
numpy.ndarray of shape (4, 4)
ddg.math.euclidean.reflection_in_a_hyperplane(normal, level)[source]

Calculate the reflection matrix (in homogeneous coordinates) of a hyperplane.

Parameters:
normalarray_like of shape (n,)

Normal of the plane

levelfloat

Level of the plane

Returns:
reflection_matrixnp.ndarray of shape (n + 1, n + 1)
Raises:
ValueError

If the given normal does not have shape (n,).

Examples

>>> import numpy as np
>>> from ddg.math.euclidean import reflection_in_a_hyperplane as re
>>> M = re((1, 0), 0)
>>> M
array([[-1.,  0.,  0.],
       [ 0.,  1.,  0.],
       [ 0.,  0.,  1.]])
>>> p = (-5, 10, 1)
>>> reflected_p = M.dot(p)
>>> reflected_p
array([ 5., 10.,  1.])
>>> M = re((1, 1, 1), 5)
>>> M
array([[ 0.33333333, -0.66666667, -0.66666667,  5.77350269],
       [-0.66666667,  0.33333333, -0.66666667,  5.77350269],
       [-0.66666667, -0.66666667,  0.33333333,  5.77350269],
       [ 0.        ,  0.        ,  0.        ,  1.        ]])
>>> p = (0, 0, 0, 1)
>>> reflected_p = M.dot(p)
>>> reflected_p
array([5.77350269, 5.77350269, 5.77350269, 1.        ])
ddg.math.euclidean.catmull_rom_spline(P0, P1, P2, P3, alpha=0.5, nPoints=100)[source]

Computes the points of the Catmull-Rom spline segment, using four control points.

Parameters:
P0array_like of shape (2,) or (3,)

First control point, which is used to determine curvature of spline segment between P1 and P2.

P1array_like of shape (2,) or (3,)

Second control point, which forms one end of the spline segment.

P2array_like of shape (2,) or (3,)

Third control point, which forms the other end of the spline segment.

P3array_like of shape (2,) or (3,)

Fourth control point, which is used to determine curvature of spline segment between P1 and P2.

alphafloat (default=0.5)

Spline Knot parameter, ranges from 0 to 1.

nPointsint (default=100)

Number of points making up the resulting spline segment.

Returns:
numpy.ndarray of shape (nPoints, 2) or (nPoints, 3)

The Catmull-Rom spline segment between points P1 and P2.

Notes

Depending on the value of alpha, we get uniform (alpha=0), centripetal (alpha=0.5) and chordal (alpha=1) parameterization of Catmull–Rom spline.

ddg.math.euclidean.catmull_rom_curve(P, alpha=0.5, nPoints=4)[source]

Calculate Catmull Rom for a chain of points and return the combined curve.

Parameters:
Parray_like of array_likes each of shape (2,) or (3,)

Chain of points from which the quadruples for the catmull_rom_spline function are taken.

alphafloat (default=0.5)

Spline Knot parameter, ranges from 0 to 1.

nPointsint (default=4)

The number of points to include in each curve segment.

Returns:
numpy.ndarray of shape (m, 2) or (m, 3)

Catmull-Rom curve made up of n-1 spline segments, each containing nPoints, i.e m = 2 - n + nPoints*(n-1) where n denotes the length of given list of points.

Notes

The resulting curve passes through all the given points.

The centripetal parameterization (alpha=0.5) of Catmull-Rom curve is the only parameterization that guarantees that the curves do not form cusps or self-intersections within its segments.

ddg.math.euclidean.extend_to_onb(vectors, index=0)[source]

Extends a given list of vectors to an orthonormal basis for n-dimensional space.

Parameters:
vectorslist of numpy.ndarray of shape (n,) or a numpy.ndarray of shape (k, n)

List of linearly independent vectors to complete into an orthonormal basis.

indexint (default=0)

If only one vector is given, it’s possible to specify its index in the resulting orthonormal basis. If more than one vector is given, the first k columns of the result will have the same span as the given family of vectors.

Returns:
numpy.ndarray of shape (n, n)

Columns of resulting matrix form an orthonormal basis.

Raises:
ValueError

If vectors are given as numpy.ndarray of dimension greater than 2.

ValueError

If the given family of vectors is linearly dependent.

Notes

The given vector(s) will be normalized.

This function uses Gram-Schmidt, so it might be a bit unstable.

If the vectors are given as a 2-dimensional numpy.ndarray, the rows of the matrix will be considered as input.

ddg.math.euclidean.distance(p0, p1)[source]

Computes the Euclidean distance between two given point in n-dimensional space.

Parameters:
p0, p1array_like of shape (n,)

Points to calculate distance of.

Returns:
float

Euclidean distance, i.e l2 norm of difference between the given two point.

ddg.math.euclidean.distance_lines(a, v, b, w)[source]

Calculate the shortest distance between the lines a + sv and b + tw. This only works in 3D-space.

Parameters:
a, bnumpy.ndarray of shape (3,)

Initial vector of line.

v, warray_like of shape (3,)

Direction vector of line.

Returns:
float

Shortest distance of line.

ddg.math.euclidean.skew_symmetric_matrix(v)[source]

Returns the skew-symmetric matrix

 0  | -v3  |  v2
----+------+----
 v3 |  0   | -v1
----+------+----
-v2 |  v1  |  0
Parameters:
vnumpy.ndarray of shape (3,)
Returns:
numpy.ndarray of shape (3, 3)
ddg.math.euclidean.normalize(v)[source]

Normalizes a given vector in n-dimensional space.

Parameters:
varray_like of shape (n,)

The vector to normalize.

Returns:
numpy.ndarray of shape (n,)

The normalized vector if it is non-zero, the zero vector otherwise.

ddg.math.euclidean.embed(v, level=0, index=-1)[source]

Embeds a coordinate vector by inserting the given level.

Parameters:
varray_like of shape (n,)

The vector to be embedded.

levelfloat (default=0)

The value to insert.

indexint (default=-1)

Index of inserted value in the output.

Returns:
numpy.ndarray of shape (n+1,)
ddg.math.euclidean.angle_bisector_orientation_preserving(n1, d1, n2, d2)[source]

Calculate the orientation preserving angle bisector of two hyperplanes.

The hyperplanes are represented by a n-dimensional normal vector and a level.

Parameters:
n1array_like of shape (n,)

Normal of first hyperplane.

d1float

Level of first hyperplane (distace to origin).

n2array_like of shape (n,)

Normal of second hyperplane.

d2float

Level of second hyperplane (distace to origin).

Returns:
n, dnumpy.ndarray of shape (n,), float

Normal and level of orientation preserving angle bisector.

ddg.math.euclidean.angle_bisector_orientation_reversing(n1, d1, n2, d2)[source]

Calculate the orientation reversing angle bisector of two hyperplanes.

The hyperplanes are represented by a n-dimensional normal vector and a level.

Parameters:
n1array_like of shape (n,)

Normal of first hyperplane.

d1float

Level of first hyperplane (distace to origin).

n2array_like of shape (n,)

Normal of second hyperplane.

d2float

Level of second hyperplane (distace to origin).

Returns:
n, dnumpy.ndarray of shape (n,), float

Normal and level of orientation reversing angle bisector.