Knot

Class for dealing with curves as knots. Knot provides many methods for topological manipulation and calculations.

API documentation

class pyknotid.spacecurves.knot.Knot(points, verbose=True, add_closure=False, zero_centroid=False)[source]

Bases: pyknotid.spacecurves.spacecurve.SpaceCurve

Class for holding the vertices of a single line, providing helper methods for convenient manipulation and analysis.

A Knot just represents a single space curve, it may be topologically trivial!

This class deliberately combines methods to do many different kinds of measurements or manipulations. Some of these are externally available through other modules in pyknotid - if so, this is usually indicated in the method docstrings.

Parameters:
  • points (array-like) – The 3d points (vertices) of a piecewise linear curve representation
  • verbose (bool) – Indicates whether the Knot should print information during processing
  • add_closure (bool) – If True, adds a final point to the knot near to the start point, so that it will appear visually to close when plotted.
alexander_at_root(root, round=True, **kwargs)[source]

Returns the Alexander polynomial at the given root of unity, i.e. evaluated at exp(2 pi I / root).

The result returned is the absolute value.

Parameters:
  • root (int) – The root of unity to use, i.e. evaluating at exp(2 pi I / root). If this is iterable, this method returns a list of the results at every value of that iterable.
  • round (bool) – If True and n in (1, 2, 3, 4), the result will be rounded to the nearest integer for convenience, and returned as an integer type.
  • **kwargs – These are passed directly to alexander_polynomial().
alexander_polynomial(variable=-1, quadrant='lr', mode='python', **kwargs)[source]

Returns the Alexander polynomial at the given point, as calculated by pyknotid.invariants.alexander().

See pyknotid.invariants.alexander() for the meanings of the named arguments.

copy()[source]

Returns another knot with the same points and verbosity as self. Other attributes (e.g. cached crossings) are not preserved.

determinant()[source]

Returns the determinant of the knot. This is the Alexander polynomial evaluated at -1.

exterior_manifold()[source]

The knot complement manifold of self as a SnapPy class giving access to all of SnapPy’s tools.

This method requires that Spherogram, and possibly SnapPy, are installed.

hyperbolic_volume()[source]

Returns the hyperbolic volume at the given point, via pyknotid.representations.PlanarDiagram.as_spherogram().

Returns:
  • volume (float) – A float representing the volume returned.
  • accuracy (int) – The number of digits of precision. This is significant digits, e.g. 0.00021 with 1 digit precision = 2E-4.
  • solution_type (str) – The solution type of the manifold. Normally one of: - ‘contains degenerate tetrahedra’ => may not be a valid result - ‘all tetrahedra positively oriented’ =>
    really probably hyperbolic
identify(determinant=True, alexander=False, roots=(2, 3, 4), min_crossings=True)[source]

Provides a simple interface to pyknotid.catalogue.identify.from_invariants(), by passing the given invariants. This does not support all invariants available, or more sophisticated identification methods, so don’t be afraid to use the catalogue functions directly.

Parameters:
  • determinant (bool) – If True, uses the knot determinant in the identification. Defaults to True.
  • alexander (bool) – If True-like, uses the full alexander polynomial in the identification. If the input is a dictionary of kwargs, these are passed straight to self.alexander_polynomial.
  • roots (iterable) – A list of roots of unity at which to evaluate. Defaults to (2, 3, 4), the first of which is redundant with the determinant. Note that higher roots can be calculated, but aren’t available in the database.
  • min_crossings (bool) – If True, the output is restricted to knots with fewer crossings than the current projection of this one. Defaults to True. The only reason to turn this off is to see what other knots have the same invariants, it is never not useful for direct identification.
isolate_knot()[source]

Return indices of self.points within which the knot (if any) appears to lie, according to a simple closure algorithm.

This method is experimental and may not provide very good results.

planar_writhe_quantities(num_angles=100, **kwargs)[source]

Returns the second order writhes, and arnold 2St+J+ values, for a range of different projection directions.

plot(**kwargs)[source]

Plots the line. See pyknotid.visualise.plot_line() for full documentation.

plot_isolated(**kwargs)[source]

Plots the curve in red, except for the isolated local knot which is coloured blue. The local knot is found with self.isolate_knot, which may not be reliable or have good resolution.

Parameters:**kwargs – kwargs are passed directly to Knot.plot().
points

The points of the spacecurve, as an Nx3 numpy array.

slipknot_alexander(num_samples=0, **kwargs)[source]
Parameters:
  • num_samples (int) – The number of indices to cut at. Defaults to 0, which means to sample at all indices.
  • **kwargs – Keyword arguments, passed directly to :meth:`pyknotid.spacecurves.openknot.OpenKnot.alexander_fractions.
vassiliev_degree_2(simplify=True, **kwargs)[source]

Returns the Vassiliev invariant of degree 2 for the Knot.

Parameters:
  • simplify (bool) – If True, simplifies the Gauss code of self before calculating the invariant. Defaults to True, but will work fine if you set it to False (and might even be faster).
  • **kwargs – These are passed directly to gauss_code().
vassiliev_degree_3(simplify=True, try_cython=True, **kwargs)[source]

Returns the Vassiliev invariant of degree 3 for the Knot.

Parameters:
  • simplify (bool) – If True, simplifies the Gauss code of self before calculating the invariant. Defaults to True, but will work fine if you set it to False (and might even be faster).
  • try_cython (bool) – Whether to try and use an optimised cython version of the routine (takes about 1/3 of the time for complex representations). Defaults to True, but the python fallback will be slower than setting it to False if the cython function is not available.
  • **kwargs – These are passed directly to gauss_code().
whitney_index()[source]

The degree of the Gauss map mapping a point on the curve to the direction of the positive tangent vector at this point.