euclidean

A loose collection of tools for Euclidean geometry in the plane.

See also

flatsurf.geometry.circle for everything specific to circles in the plane

flatsurf.geometry.euclidean.acos(cos_angle, numerical=False)[source]

Return the arccosine of cos_angle as a multiple of 2π, i.e., as a value between 0 and 1/2.

INPUT:

  • cos_angle – a floating point number, the cosine of an angle

  • numerical – a boolean (default: False); whether to return a numerical approximation of the arccosine or try to reconstruct an exact rational value for the arccosine (in radians).

EXAMPLES:

sage: from flatsurf.geometry.euclidean import acos

sage: acos(1)
0
sage: acos(.5)
1/6
sage: acos(0)
1/4
sage: acos(-.5)
1/3
sage: acos(-1)
1/2

sage: acos(.25)
Traceback (most recent call last):
...
NotImplementedError: cannot recover a rational angle from these numerical results
sage: acos(.25, numerical=True)
0.2097846883724169
flatsurf.geometry.euclidean.angle(u, v, numerical=False)[source]

Return the angle between the vectors u and v divided by \(2 \pi\).

INPUT:

  • u, v - vectors

  • numerical - boolean (default: False), whether to return floating point numbers

EXAMPLES:

sage: from flatsurf.geometry.euclidean import angle

As the implementation is dirty, we at least check that it works for all denominator up to 20:

sage: u = vector((AA(1),AA(0)))
sage: for n in xsrange(1,20):       # long time  (1.5s)
....:     for k in xsrange(1,n):
....:         v = vector((AA(cos(2*k*pi/n)), AA(sin(2*k*pi/n))))
....:         assert angle(u,v) == k/n

The numerical version (working over floating point numbers):

sage: import math
sage: u = (1, 0)
sage: for n in xsrange(1,20):
....:     for k in xsrange(1,n):
....:         a = 2 * k * math.pi / n
....:         v = (math.cos(a), math.sin(a))
....:         assert abs(angle(u,v,numerical=True) * 2 * math.pi - a) < 1.e-10

If the angle is not rational, then the method returns an element in the real lazy field:

sage: v = vector((AA(sqrt(2)), AA(sqrt(3))))
sage: a = angle(u, v)
Traceback (most recent call last):
...
NotImplementedError: cannot recover a rational angle from these numerical results
sage: a = angle(u, v, numerical=True)
sage: a    # abs tol 1e-14
0.14102355421224375
sage: exp(2*pi.n()*CC(0,1)*a)
0.632455532033676 + 0.774596669241483*I
sage: v / v.norm()
(0.6324555320336758?, 0.774596669241484?)
flatsurf.geometry.euclidean.ccw(v, w)[source]

Return a positive number if the turn from v to w is counterclockwise, a negative number if it is clockwise, and zero if the two vectors are collinear.

Note

This function is sometimes also referred to as the wedge product or simply the determinant. We chose the more customary name ccw from computational geometry here.

EXAMPLES:

sage: from flatsurf.geometry.euclidean import ccw
sage: ccw((1, 0), (0, 1))
1
sage: ccw((1, 0), (-1, 0))
0
sage: ccw((1, 0), (0, -1))
-1
sage: ccw((1, 0), (1, 0))
0
flatsurf.geometry.euclidean.is_anti_parallel(v, w)[source]

Return whether the vectors v and w are anti-parallel, i.e., whether v and -w are parallel.

EXAMPLES:

sage: from flatsurf.geometry.euclidean import is_anti_parallel
sage: V = QQ**2

sage: is_anti_parallel(V((0,1)), V((0,-2)))
True
sage: is_anti_parallel(V((1,-1)), V((-2,2)))
True
sage: is_anti_parallel(V((4,-2)), V((-2,1)))
True
sage: is_anti_parallel(V((-1,-2)), V((2,4)))
True

sage: is_anti_parallel(V((1,1)), V((1,2)))
False
sage: is_anti_parallel(V((1,2)), V((2,1)))
False
sage: is_anti_parallel(V((0,2)), V((0,1)))
False
sage: is_anti_parallel(V((1,2)), V((1,-2)))
False
sage: is_anti_parallel(V((1,2)), V((-1,2)))
False
sage: is_anti_parallel(V((2,-1)), V((-2,-1)))
False
flatsurf.geometry.euclidean.is_between(e0, e1, f)[source]

Check whether the vector f is strictly in the sector formed by the vectors e0 and e1 (in counter-clockwise order).

EXAMPLES:

sage: from flatsurf.geometry.euclidean import is_between
sage: is_between((1, 0), (1, 1), (2, 1))
True

sage: from itertools import product
sage: vecs = [(1, 0), (1, 1), (0, 1), (-1, 1), (-1, 0), (-1, -1), (0, -1), (1, -1)]
sage: for (i, vi), (j, vj), (k, vk) in product(enumerate(vecs), repeat=3):
....:     assert is_between(vi, vj, vk) == ((i == j and i != k) or i < k < j or k < j < i or j < i < k), ((i, vi), (j, vj), (k, vk))
flatsurf.geometry.euclidean.is_cosine_sine_of_rational(cos, sin, scaled=False)[source]

Check whether the given pair is a cosine and sine of a same rational angle.

INPUT:

  • cos – a number

  • sin – a number

  • scaled – a boolean (default: False); whether to allow cos and sin to be scaled by the same positive algebraic number

EXAMPLES:

sage: from flatsurf.geometry.euclidean import is_cosine_sine_of_rational

sage: c = s = AA(sqrt(2))/2
sage: is_cosine_sine_of_rational(c, s)
True

sage: c = AA(sqrt(3))/2
sage: s = AA(1/2)
sage: is_cosine_sine_of_rational(c, s)
True

sage: c = AA(sqrt(5)/2)
sage: s = (1 - c**2).sqrt()
sage: c**2 + s**2
1.000000000000000?
sage: is_cosine_sine_of_rational(c, s)
False

sage: c = (AA(sqrt(5)) + 1)/4
sage: s = (1 - c**2).sqrt()
sage: is_cosine_sine_of_rational(c, s)
True

sage: K.<sqrt2> = NumberField(x**2 - 2, embedding=1.414)
sage: is_cosine_sine_of_rational(K.zero(), -K.one())
True
flatsurf.geometry.euclidean.is_parallel(v, w)[source]

Return whether the vectors v and w are parallel (but not anti-parallel).

EXAMPLES:

sage: from flatsurf.geometry.euclidean import is_parallel
sage: is_parallel((0, 1), (0, 1))
True
sage: is_parallel((0, 1), (0, 2))
True
sage: is_parallel((0, 1), (0, -2))
False
sage: is_parallel((0, 1), (0, 0))
False
sage: is_parallel((0, 1), (1, 0))
False
flatsurf.geometry.euclidean.is_segment_intersecting(e1, e2)[source]

Return whether the segments e1 and e2 intersect.

OUTPUT:

  • 0 - do not intersect

  • 1 - one endpoint in common

  • 2 - non-trivial intersection

EXAMPLES:

sage: from flatsurf.geometry.euclidean import is_segment_intersecting
sage: is_segment_intersecting(((0,0),(1,0)),((0,1),(0,3)))
0
sage: is_segment_intersecting(((0,0),(1,0)),((0,0),(0,3)))
1
sage: is_segment_intersecting(((0,0),(1,0)),((0,-1),(0,3)))
2
sage: is_segment_intersecting(((-1,-1),(1,1)),((0,0),(2,2)))
2
sage: is_segment_intersecting(((-1,-1),(1,1)),((1,1),(2,2)))
1
flatsurf.geometry.euclidean.line_intersection(p1, p2, q1, q2)[source]

Return the point of intersection between the line joining p1 to p2 and the line joining q1 to q2. If the lines do not have a single point of intersection, we return None. Here p1, p2, q1 and q2 should be vectors in the plane.

flatsurf.geometry.euclidean.projectivization(x, y, signed=True, denominator=None)[source]

Return a simplified version of the projective coordinate [x: y].

If signed (the default), the second coordinate is made non-negative; otherwise the coordinates keep their signs.

If denominator is False, returns [x/y: 1] up to sign. Otherwise, the returned coordinates have no denominator and no non-unit gcd.

flatsurf.geometry.euclidean.slope(a, rotate=1)[source]

Return either 1 (positive slope) or -1 (negative slope).

If rotate is set to 1 then consider the edge as if it was rotated counterclockwise infinitesimally.

EXAMPLES:

sage: from flatsurf.geometry.euclidean import slope
sage: slope((1, 1))
1
sage: slope((-1, 1))
-1
sage: slope((-1, -1))
1
sage: slope((1, -1))
-1

sage: slope((1, 0))
1
sage: slope((0, 1))
-1
sage: slope((-1, 0))
1
sage: slope((0, -1))
-1

sage: slope((1, 0), rotate=-1)
-1
sage: slope((0, 1), rotate=-1)
1
sage: slope((-1, 0), rotate=-1)
-1
sage: slope((0, -1), rotate=-1)
1

sage: slope((1, 0), rotate=0)
0
sage: slope((0, 1), rotate=0)
0
sage: slope((-1, 0), rotate=0)
0
sage: slope((0, -1), rotate=0)
0

sage: slope((0, 0))
Traceback (most recent call last):
...
ValueError: zero vector
flatsurf.geometry.euclidean.solve(x, u, y, v)[source]

Return (a,b) so that: x + au = y + bv

INPUT:

  • x, u, y, v – two dimensional vectors

EXAMPLES:

sage: from flatsurf.geometry.euclidean import solve
sage: K.<sqrt2> = NumberField(x^2 - 2, embedding=AA(2).sqrt())
sage: V = VectorSpace(K,2)
sage: x = V((1,-sqrt2))
sage: y = V((1,1))
sage: a = V((0,1))
sage: b = V((-sqrt2, sqrt2+1))
sage: u = V((0,1))
sage: v = V((-sqrt2, sqrt2+1))
sage: a, b = solve(x,u,y,v)
sage: x + a*u == y + b*v
True

sage: u = V((1,1))
sage: v = V((1,sqrt2))
sage: a, b = solve(x,u,y,v)
sage: x + a*u == y + b*v
True