SpaceCurve¶
The SpaceCurve
class is the base for all space-curve analysis
in pyknotid. It provides methods for geometrical manipulation
(translation, rotation etc.), calculating geometrical characteristics
such as the writhe, and obtaining topological representations of the
curve by analysing its crossings in projection.
pyknotid provides other classes for topological analysis of the curves:
Knot
for calculating knot invariants of the space-curve.Link
to handle multiple SpaceCurves and calculate linking invariants.Cell
for handling multiple space-curves in a box with periodic boundaries.
API documentation¶
-
class
pyknotid.spacecurves.spacecurve.
SpaceCurve
(points, verbose=True, add_closure=False, zero_centroid=False)[source]¶ Bases:
object
Class for holding the vertices of a single line, providing helper methods for convenient manipulation and analysis.
The methods of this class are largely geometrical (though this includes listing the crossings in projection and extracting a Gauss code etc.). For topological measurements, you should use a
Knot
.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 SpaceCurve 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.
- zero_centroid (bool) – If True, shifts the coordinates of the points so that their centre of mass is at the origin.
-
arclength
(include_closure=True)[source]¶ Returns the arclength of the line, the sum of lengths of each piecewise linear segment.
Parameters: include_closure (bool) – Whether to include the distance between the final and first points. Defaults to True.
-
average_crossing_number
(samples=10, recalculate=False, **kwargs)[source]¶ The (approximate) average crossing number of the space curve, obtained by averaging the planar writhe over the given number of directions.
Parameters: - samples (int) – The number of directions to average over.
- recalculate (bool) – Whether to recalculate the ACN.
- **kwargs – These are passed directly to
raw_crossings()
.
-
close
()[source]¶ Adds the starting point to the end of the curve, so that it ends exactly where it began.
-
classmethod
closing_on_sphere
(line, com=(0.0, 0.0, 0.0))[source]¶ Adds new vertices to close the line at its maximum radius, returning a SpaceCurve representing the result.
Parameters: - line (ndarray) – The points of the line.
- com (iterable) – Optional additional centre of mass to shift by before closing
-
copy
()[source]¶ Returns another knot with the same points and verbosity as self. Other attributes (e.g. cached crossings) are not preserved.
-
cuaps
(include_closure=True)[source]¶ Returns a list of the ‘cuaps’, where the curve is parallel to the positive x-axis. See D Bar-Natan and R van der Veen, “A polynomial time knot polynomial”, 2017.
-
classmethod
from_braid_word
(word)[source]¶ Returns a
SpaceCurve
instance formed from the given braid word.The braid word should be of the form ‘aAbBcC’ (i.e. capitalisation denotes inverse).
Parameters: word (str) – The braid word to interpret.
-
classmethod
from_csv
(filen, **kwargs)[source]¶ Loads knot points from the given csv file, and returns a
SpaceCurve
with those points.Arguments are passed straight to
pyknot.io.from_csv()
.
-
classmethod
from_gauss_code
(code)[source]¶ Creates a Knot from the given code, which must be provided as a string and may optionally include crossing orientations (these are actually ignored).
-
classmethod
from_json
(filen)[source]¶ Loads knot points from the given filename, assuming json format, and returns a
SpaceCurve
with those points.
-
classmethod
from_lattice_data
(line)[source]¶ Returns a
SpaceCurve
instance in which the line has been slightly translated and rotated, in order to (practically) ensure no self intersections in closure or coincident points in projection.Parameters: line (array-like) – The list of points in the line. May be any type that SpaceCurve normally accepts. Returns: Return type: SpaceCurve
-
classmethod
from_periodic_line
(line, shape, perturb=True, **kwargs)[source]¶ Returns a
SpaceCurve
instance in which the line has been unwrapped through the periodic boundaries.Parameters: - line (array-like) – The Nx3 vector of points in the line
- shape (array-like) – The x, y, z distances of the periodic boundary
- perturb (bool) – If True, translates and rotates the knot to avoid any lattice problems.
-
gauss_code
(recalculate=False, **kwargs)[source]¶ Returns a
GaussCode
instance representing the crossings of the knot.The GaussCode instance is cached internally. If you want to recalculate it (e.g. to get an unsimplified version if you have simplified it), you should pass recalculate=True.
This method passes kwargs directly to
raw_crossings()
, see the documentation of that function for all options.
-
gauss_diagram
(simplify=False, **kwargs)[source]¶ Returns a
GaussDiagram
instance representing the crossings of the knot.This method passes kwargs directly to
raw_crossings()
, see the documentation of that function for all options.
-
interpolate
(num_points, s=0, **kwargs)[source]¶ Replaces self.points with points from a B-spline interpolation.
This method uses scipy.interpolate.splprep. kwargs are passed to this function.
-
octree_simplify
(runs=1, plot=False, rotate=True, obey_knotting=True, **kwargs)[source]¶ Simplifies the curve via the octree reduction of :module:`pyknotid.simplify.octree`.
Parameters: - runs (int) – The number of times to run the octree simplification. Defaults to 1.
- plot (bool) – Whether to plot the curve after each run. Defaults to False.
- rotate (bool) – Whether to rotate the space curve before each run. Defaults to True as this can make things much faster.
- obey_knotting (bool) – Whether to not let the line pass through itself. Defaults to True as this is always what you want for a closed curve.
- **kwargs – Any remaining kwargs are passed to the
pyknotid.simplify.octree.OctreeCell
constructor.
-
planar_diagram
(**kwargs)[source]¶ Returns a
PlanarDiagram
instance representing the crossings of the knot.This method passes kwargs directly to
raw_crossings()
, see the documentation of that function for all options.
-
planar_second_order_writhe
(**kwargs)[source]¶ The second order writhe (type 2, i1,i3,i2,i4) of the projection of the curve along the z axis.
-
planar_writhe
(**kwargs)[source]¶ Returns the current planar writhe of the knot; the signed sum of crossings of the current projection.
The ‘true’ writhe is the average of this quantity over all projection directions, and is available from the
writhe()
method.Parameters: **kwargs – These are passed directly to raw_crossings()
.
-
plot
(mode='auto', clf=True, closed=False, **kwargs)[source]¶ Plots the line. See
pyknotid.visualise.plot_line()
for full documentation.
-
plot_projection
(with_crossings=True, mark_start=False, fig_ax=None, show=True, mark_points=False)[source]¶ Plots a 2D diagram of the knot projected along the current z-axis. The crossings, and start point of the curve, can optionally be marked.
The projection is drawn using matplotlib.
Parameters: - with_crossings (bool) – If True, marks the location of each self-intersection in projection. Defaults to True.
- mark_start (bool) – If True, marks the first point of the curve. Default to False.
- fig_ax (tuple) – A 2-tuple of the matplotlib (fig, ax) to use, or None to create a new pair.
- show (bool) – If True, opens a new window showing the drawing. Defaults to True.
Returns: Return type: A 2-tuple of the matplotlib (fig, ax) used for the drawing.
-
points
¶ The points of the spacecurve, as an Nx3 numpy array.
-
radius_of_gyration
()[source]¶ Returns the radius of gyration of the points of self, assuming each has equal weight and ignoring the connecting lines.
-
raw_crossings
(mode='use_max_jump', include_closure=True, recalculate=False, try_cython=True)[source]¶ Returns the crossings in the diagram of the projection of the space curve into its z=0 plane.
The crossings will be calculated the first time this function is called, then cached until an operation that would change the list (e.g. rotation, or changing
self.points
).Multiple modes are available (see parameters) - you should be aware of this because different modes may be vastly slower or faster depending on the type of line.
Parameters: - mode (str, optional) – One of
'count_every_jump'
or'use_max_jump'
or'naive'
. In the first case, walking along the line uses information about the length of every step. In the second, it guesses that all steps have the same length as the maximum step length. In the last, no optimisation is made and every crossing is checked. The optimal choice depends on the data but is usually'use_max_jump'
, which is the default. - include_closure (bool, optional) – Whether to include crossings with the line joining the start and end points. Defaults to True.
- recalculate (bool, optional) – Whether to force a recalculation of the crossing positions. Defaults to False.
- try_cython (bool, optional) – Whether to force the use of the python (as opposed to cython) implementation of find_crossings. This will make no difference if the cython could not be loaded, in which case python is already used automatically. Defaults to True.
Returns: The raw array of floats representing crossings, of the form [[line_index, other_index, +-1, +-1], …], where the line_index and other_index are in arclength parameterised by integers for each vertex and linearly interpolated, and the +-1 represent over/under and clockwise/anticlockwise respectively.
Return type: array-like
- mode (str, optional) – One of
-
reparameterised
(mode='arclength', num_points=None, interpolation='linear')[source]¶ Returns a new
SpaceCurve
where new points have been selected by interpolating the current ones.Warning
This doesn’t do what you expect! The new segments will probably not all be separated by the right amount in terms of the new parameterisation.
Parameters: - mode (str) – The function to reparameterise by. Defaults to ‘arclength’, which is currently the only option.
- num_points (int) – The number of points in the new parameterisation. Defaults to None, which means the same as the current number.
- interpolation (str) – The type of interpolation to use, passed directly to the
kind
option ofscipy.interpolate.interp1d
. Defaults to ‘linear’, and other options have not been tested.
-
representation
(recalculate=False, **kwargs)[source]¶ Returns a
Representation
instance representing the crossings of the knot.The Representation instance is cached internally. If you want to recalculate it (e.g. to get an unsimplified version if you have simplified it), you should pass recalculate=True.
This method passes kwargs directly to
raw_crossings()
, see the documentation of that function for all options.
-
rotate
(angles=None)[source]¶ Rotates all the points of self by the given angles in each axis.
Parameters: angles (array-like) – The rotation angles about x, y and z. If None, random angles are used. Defaults to None.
-
scale
(factor)[source]¶ Scales all the points of self by the given factor.
You can accomplish the same thing, or other more subtle transformations, by modifying self.:py:attr:points.
-
segment_arclengths
()[source]¶ Returns an array of arclengths of every step in the line defined by self.points.
-
simplify_straight_segments
(closed=False)[source]¶ Replaces successive curve segments with identical tangents with a single longer segment.
-
smooth
(repeats=1, periodic=True, window_len=10, window='hanning')[source]¶ Smooths each of the x, y and z components of self.points by convolving with a window of the given type and size.
Warning
This is not topologically safe, it can change the knot type of the curve. For topologically safe reduction, see
octree_simplify()
.Parameters: - repeats (int) – Number of times to run the smoothing algorithm. Defaults to 1.
- periodic (bool) – If True, the convolution window wraps around the curve. Defaults to True.
- window_len (int) – Width of the convolution window. Defaults to 10.
Passed to
pyknotid.spacecurves.smooth.smooth()
. - window (string) – The type of convolution window. Defaults to ‘hanning’.
Passed to
pyknotid.spacecurves.smooth.smooth()
.
-
to_json
(filen)[source]¶ Writes the knot points to the given filename, in a json format that can be read later by
SpaceCurve.from_json()
. Usespyknotid.io.to_json_file()
internally.
-
to_txt
(filen)[source]¶ Writes the knot points to the given filename, formatted with each x,y,z component of each point space-separated on its own line, i.e.:
... 1.2 6.1 98.5 6.19 8.5 1.9 ...
-
translate
(vector)[source]¶ Translates all the points of self by the given vector.
Parameters: vector (array-like) – The x, y, z translation distances
-
writhe
(samples=10, recalculate=False, method='integral', include_acn=False, **kwargs)[source]¶ The (approximate) writhe of the space curve, obtained by averaging the planar writhe over the given number of directions.
Parameters: - samples (int) – The number of directions to average over. Defaults to 10.
- recalculate (bool) – Whether to recalculate the writhe even if a cached result is available. Defaults to False.
- method (str) – If ‘projections’, averages the planar writhe over many projections. If ‘integral’, calculates the writhing integral.
- **kwargs – These are passed directly to
raw_crossings()
.