finitely_generated_matrix_group

Class for matrix groups generated by a finite number of elements.

EXAMPLES:

sage: from flatsurf.geometry.finitely_generated_matrix_group import  FinitelyGenerated2x2MatrixGroup

sage: m1 = matrix([[1,1],[0,1]])
sage: m2 = matrix([[1,0],[1,1]])
sage: G = FinitelyGenerated2x2MatrixGroup([m1,m2])
sage: G
Matrix group generated by:
[1 1]  [1 0]
[0 1], [1 1]
sage: it = iter(G)
sage: [next(it) for _ in range(5)]
[
[1 0]  [1 1]  [1 2]  [2 1]  [ 0  1]
[0 1], [0 1], [0 1], [1 1], [-1  1]
]

sage: G = FinitelyGenerated2x2MatrixGroup([identity_matrix(2)])
class flatsurf.geometry.finitely_generated_matrix_group.FinitelyGenerated2x2MatrixGroup(matrices, matrix_space=None, category=None)[source]

Finitely generated group of 2x2 matrices with real coefficients

See also

sage.groups.group for the general interface of groups like this in SageMath

cardinality()[source]
gen(i)[source]
is_abelian()[source]

Check whether this group is abelian.

EXAMPLES:

sage: from flatsurf.geometry.finitely_generated_matrix_group import  FinitelyGenerated2x2MatrixGroup

sage: m1 = matrix([[1,1],[0,1]])
sage: m2 = matrix([[1,0],[1,1]])
sage: G = FinitelyGenerated2x2MatrixGroup([m1,m2])
sage: G.is_abelian()
False

sage: G = FinitelyGenerated2x2MatrixGroup([m1])
sage: G.is_abelian()
True
is_finite()[source]

Check whether the group is finite.

A group is finite if and only if it is conjugate to a (finite) subgroup of O(2). This is actually also true in higher dimensions.

EXAMPLES:

sage: from flatsurf.geometry.finitely_generated_matrix_group import FinitelyGenerated2x2MatrixGroup
sage: G = FinitelyGenerated2x2MatrixGroup([identity_matrix(2)])
sage: G.is_finite()
True

sage: t = matrix(2, [2,1,1,1])

sage: m1 = matrix([[0,1],[-1,0]])
sage: m2 = matrix([[1,-1],[1,0]])
sage: FinitelyGenerated2x2MatrixGroup([m1]).is_finite()
True
sage: FinitelyGenerated2x2MatrixGroup([t*m1*~t]).is_finite()
True
sage: FinitelyGenerated2x2MatrixGroup([m2]).is_finite()
True
sage: FinitelyGenerated2x2MatrixGroup([m1,m2]).is_finite()
False
sage: FinitelyGenerated2x2MatrixGroup([t*m1*~t,t*m2*~t]).is_finite()
False

sage: from flatsurf.geometry.polygon import number_field_elements_from_algebraics
sage: c5 = QQbar.zeta(5).real()
sage: s5 = QQbar.zeta(5).imag()
sage: K, (c5,s5) = number_field_elements_from_algebraics([c5,s5])
sage: r = matrix(K, 2, [c5,-s5,s5,c5])
sage: FinitelyGenerated2x2MatrixGroup([m1,r]).is_finite()
True
sage: FinitelyGenerated2x2MatrixGroup([t*m1*~t,t*r*~t]).is_finite()
True
sage: FinitelyGenerated2x2MatrixGroup([m2,r]).is_finite()
False
sage: FinitelyGenerated2x2MatrixGroup([t*m2*~t, t*r*~t]).is_finite()
False
one()[source]
some_elements()[source]
flatsurf.geometry.finitely_generated_matrix_group.contains_definite_form(V)[source]

Check whether a given a subspace of the 3 dimensional space (a,b,c) contains a definitive positive quadratic form ax^2 + 2bxy + by^2.

flatsurf.geometry.finitely_generated_matrix_group.invariant_quadratic_forms(m)[source]

Return the space of quadratic forms invariant under m.

If the matrix is orientable, there is only one (what about strange eigenvalues?)

If the det(m) == -1 and trace(m) == 0, there are 2.

EXAMPLES:

sage: from flatsurf.geometry.finitely_generated_matrix_group import invariant_quadratic_forms

sage: invariant_quadratic_forms(matrix(2, [0,1,-1,0]))
Free module of degree 3 and rank 1 over Integer Ring
Echelon basis matrix:
[1 1 0]
sage: invariant_quadratic_forms(matrix(2, [1,1,0,1]))
Free module of degree 3 and rank 1 over Integer Ring
Echelon basis matrix:
[0 1 0]
sage: invariant_quadratic_forms(matrix(2, [2,1,1,1]))
Free module of degree 3 and rank 1 over Integer Ring
Echelon basis matrix:
[ 2 -2 -1]

sage: r = matrix(2,[1,-1,1,0])
sage: q = invariant_quadratic_forms(r)
sage: q
Free module of degree 3 and rank 1 over Integer Ring
Echelon basis matrix:
[ 2  2 -1]
sage: a,c,b = q.gen()
sage: m = matrix(2, [a,b,b,c])
sage: r.transpose() * m * r == m
True

sage: invariant_quadratic_forms(matrix(2, [-1,0,0,1]))
Free module of degree 3 and rank 2 over Integer Ring
Echelon basis matrix:
[1 0 0]
[0 1 0]

sage: invariant_quadratic_forms(-identity_matrix(2))
Traceback (most recent call last):
...
ValueError: m must be non scalar

sage: for _ in range(100):
....:     r = random_matrix(ZZ, 2, algorithm='unimodular')
....:     if r.is_scalar(): continue
....:     q = invariant_quadratic_forms(r)
....:     a,c,b = q.random_element()
....:     m = matrix(2, [a,b,b,c])
....:     assert r.transpose() * m * r == m
....:     m[0,0] = -m[0,0]
....:     m[0,1] = -m[0,1]
....:     q = invariant_quadratic_forms(r)
....:     a,c,b = q.random_element()
....:     m = matrix(2, [a,b,b,c])
....:     assert r.transpose() * m * r == m
flatsurf.geometry.finitely_generated_matrix_group.matrix_multiplicative_order(m)[source]

Return the order of the 2x2 matrix m.