polyhedra

Initialize self. See help(type(self)) for accurate signature.

class flatsurf.geometry.polyhedra.ConeSurfaceToPolyhedronMap(cone_surface, polyhedron, mapping_data)[source]

A map sending objects defined on a ConeSurface built from a polyhedron to the polyhedron. Currently, this works to send a trajectory to a list of points.

This class should not be called directly. You get an object of this type from polyhedron_to_cone_surface.

flatsurf.geometry.polyhedra.platonic_cube()[source]

Produce a triple consisting of a polyhedral version of the platonic cube, the associated cone surface, and a ConeSurfaceToPolyhedronMap from the surface to the polyhedron.

EXAMPLES:

sage: from flatsurf.geometry.polyhedra import platonic_cube
sage: polyhedron,surface,surface_to_polyhedron = platonic_cube()
sage: TestSuite(surface).run()
flatsurf.geometry.polyhedra.platonic_dodecahedron()[source]

Produce a triple consisting of a polyhedral version of the platonic dodecahedron, the associated cone surface, and a ConeSurfaceToPolyhedronMap from the surface to the polyhedron.

EXAMPLES:

sage: from flatsurf.geometry.polyhedra import platonic_dodecahedron
sage: polyhedron, surface, surface_to_polyhedron = platonic_dodecahedron()  # long time (1s)
sage: TestSuite(surface).run()  # long time (.8s)
flatsurf.geometry.polyhedra.platonic_icosahedron()[source]

Produce a triple consisting of a polyhedral version of the platonic icosahedron, the associated cone surface, and a ConeSurfaceToPolyhedronMap from the surface to the polyhedron.

EXAMPLES:

sage: from flatsurf.geometry.polyhedra import platonic_icosahedron
sage: polyhedron,surface,surface_to_polyhedron = platonic_icosahedron()  # long time (.9s)
sage: TestSuite(surface).run()  # long time (see above)
flatsurf.geometry.polyhedra.platonic_octahedron()[source]

Produce a triple consisting of a polyhedral version of the platonic octahedron, the associated cone surface, and a ConeSurfaceToPolyhedronMap from the surface to the polyhedron.

EXAMPLES:

sage: from flatsurf.geometry.polyhedra import platonic_octahedron
sage: polyhedron,surface,surface_to_polyhedron = platonic_octahedron()  # long time (.3s)
sage: TestSuite(surface).run()  # long time (see above)
flatsurf.geometry.polyhedra.platonic_tetrahedron()[source]

Produce a triple consisting of a polyhedral version of the platonic tetrahedron, the associated cone surface, and a ConeSurfaceToPolyhedronMap from the surface to the polyhedron.

EXAMPLES:

sage: from flatsurf.geometry.polyhedra import platonic_tetrahedron
sage: polyhedron,surface,surface_to_polyhedron = platonic_tetrahedron()
sage: TestSuite(surface).run()
flatsurf.geometry.polyhedra.polyhedron_to_cone_surface(polyhedron, use_AA=False, scaling_factor=1)[source]

Construct the Euclidean Cone Surface associated to the surface of a polyhedron and a map from the cone surface to the polyhedron.

INPUT:

  • polyhedron – A 3-dimensional polyhedron, which should be defined over something that coerces into AA

  • use_AA – If True, the surface returned will be defined over AA. If false, the algorithm will find the smallest NumberField and write the field there.

  • scaling_factor – The surface returned will have a metric scaled by multiplication by this factor (compared with the original polyhendron). This can be used to produce a surface defined over a smaller NumberField.

OUTPUT:

A pair consisting of a ConeSurface and a ConeSurfaceToPolyhedronMap.

EXAMPLES:

sage: from flatsurf.geometry.polyhedra import Polyhedron, polyhedron_to_cone_surface
sage: vertices=[]
sage: for i in range(3):
....:     temp=vector([1 if k==i else 0 for k in range(3)])
....:     for j in range(-1,3,2):
....:         vertices.append(j*temp)
sage: octahedron=Polyhedron(vertices=vertices)
sage: surface,surface_to_octahedron = \
....:     polyhedron_to_cone_surface(octahedron,scaling_factor=AA(1/sqrt(2)))
sage: TestSuite(surface).run()
sage: TestSuite(surface_to_octahedron).run()  # long time (.4s)
sage: len(surface.polygons())
8
sage: surface.base_ring()
Number Field in a with defining polynomial y^2 - 3 with a = 1.732050807568878?
sage: sqrt3=surface.base_ring().gen()
sage: tangent_bundle=surface.tangent_bundle()
sage: v=tangent_bundle(0,(0,0),(sqrt3,2))
sage: traj=v.straight_line_trajectory()
sage: traj.flow(10)
sage: traj.is_saddle_connection()
True
sage: traj.combinatorial_length()
8
sage: path3d = surface_to_octahedron(traj)
sage: len(path3d)
9
sage: # We will show that the length of the path is sqrt(42):
sage: total_length = 0
sage: for i in range(8):
....:     start = path3d[i]
....:     end = path3d[i+1]
....:     total_length += (vector(end)-vector(start)).norm()
sage: ZZ(total_length**2)
42