forgi 2.0.0 documentation

«  forgi.projection.hausdorff module   ::   Contents   ::   forgi.threedee package  »

forgi.projection.projection2d module

class forgi.projection.projection2d.Projection2D(cg, proj_direction=None, rotation=0, project_virtual_atoms=False, project_virtual_residues=[])[source]

Bases: future.types.newobject.newobject

A 2D Projection of a CoarseGrainRNA unto a 2D-plane

Parameters:
  • cg

    a CoarseGrainRNA object with 3D coordinates for every element

    Note

    The projection is generated from this cg, but it is not associated with it after construction. Thus future changes of the cg are not reflected in the projection.

  • proj_direction – a carthesian vector (in 3D space) in the direction of projection. The length of this vector is not used. If proj_direction is None, cg.project_from is used. If proj_direction and cg.project_from is None, an error is raised.
  • rotate – Degrees. Rotate the projection by this amount.
condense(cutoff)[source]

Condenses points that are within the cutoff of another point or edge (=line segment) into a new point. This function modifies this Projection2D object.

The contraction of two points works like documented in self.condense_points(self, cutoff). If a node is close to a line segment, a new point is generated between this node and the line segment. Then the original line and the original node are deleted and all connections attached to the new point.

Note

The result depends on the ordering of the dictionary holding the nodes and might thus be pseudorandomly

Parameters:cutoff – Two point with a distance smaller than cuttoff are contracted. A value below 20 is reasonable.
condense_points(cutoff=1)[source]

Condenses several projection points that are within a range of less than cutoff into one point. This function modifies this Projection2D object.

Note

The result depends on the ordering of the dictionary holding the nodes and might thus be pseudorandomly

Parameters:cutoff – Two point with a distance smaller than cuttoff are contracted. A value below 20 is reasonable.
crossingPoints

All points, where 2 segments intersect

A list of triples (key1, key2, coordinate), where coordinate is a 2D vector.

get_bounding_box(margin=0.0)[source]

Returns the coordinates for a box that contains all points of the 2D projection.

Parameters:margin – increase the bounding box in every direction by this margin.
Returns:left, right, bottom, top
get_bounding_square(margin=0.0)[source]

Returns the coordinates for a square that contains all points of the 2D projection.

Parameters:margin – increase the bounding box in every direction by this margin.
Returns:left, right, bottom, top
get_branchpoint_count(degree=None)[source]

Returns the number of branchpoint.

Note

This measure is sensitive to the resolution of a projection. In an AFM image, one might not see all branching points.

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

Parameters:degree – If degree is None, count all points with degree>=3 Else: only count (branch)points of the given degree
get_cyclebasis_len()[source]

Returns the number of cycles of length>1 in the cycle basis.

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

get_leaf_leaf_distances()[source]

Get a list of distances between any pair of leaf nodes. The distances are measured in direct line, not along the path

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

Returns:a list of floats (lengths in Angstrom)
get_longest_arm_length()[source]

Get the length of the longest arm.

An arm is a simple path from a node of degree!=2 to a node of degree 1, if all the other nodes on the path have degree 2.

Note

This measure is sensitive to the resolution of a projection the same way the length of a coastline is sensitive to the resolution.

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

Returns:The length and a tuple of points (leaf_node, corresponding_branch_point)
get_maximal_path_length()[source]

Get the maximal path length from all simple paths that traverses the projection graph from one leave node to another.

Note

This measure is sensitive to the resolution of a projection the same way the length of a coastline is sensitive to the resoltuion.

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

get_some_leaf_leaf_distances()[source]

Get a list of distances between some pairs of leaf nodes. The distances are measured in direct line, not along the path

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

Returns:a list of floats (lengths in Angstrom)
get_total_length()[source]

Returns the sum of the lengths of all edges in the projection graph.

Note

This measure is sensitive to the resolution of a projection. In an AFM image, one might not see all cycles.

Warning

Whether this code will stay in the library or not depends on future evaluation of the usefulness of this and similar descriptors.

get_vres_by_position(pos)[source]
Parameters:pos – The nucleotide position in the sequence
longest_axis = None

The longest distance between any two points of the projection.

plot(ax=None, show=False, margin=5, linewidth=None, add_labels=False, line2dproperties={}, xshift=0, yshift=0, show_distances=[], print_distances=False, virtual_atoms=True)[source]

Plots the 2D projection.

This uses modified copy-paste code by Syrtis Major (c)2014-2015 under the BSD 3-Clause license and code from matplotlib under the PSF license.

Parameters:
  • ax – The axes to draw to. You can get it by calling fig, ax=matplotlib.pyplot.subplots()
  • show – If true, the matplotlib.pyplot.show() will be called at the end of this function.
  • margin – A numeric value. The margin around the plotted projection inside the (sub-)plot.
  • linewidth – The width of the lines projection.
  • add_labels – Display the name of the corresponding coarse grain element in the middle of each segment in the projection. Either a bool or a set of labels to display.
  • line2dproperties – A dictionary. Will be passed as **kwargs to the constructor of matplotlib.lines.Line2D. See http://matplotlib.org/api/lines_api.html#matplotlib.lines.Line2D
  • yshift (xshift,) – Shift the projection by the given amount inside the canvas.
  • show_distances – A list of tuples of strings, e.g. [("h1","h8"),("h2","m15")]. Show the distances between these elements in the plot
  • print_distances – Bool. Print all distances from show_distances at the side of the plot instead of directly next to the distance
points

All points that are at the ends of coarse-grain elements. This does not include points where coarse grain elements intersect in the projection.

Some points might be duplicates

Returns:A generator yielding all points.
proj_direction

A vector describing the direction of the projection

proj_graph

A graph describing the projected RNA.

This graph is stored as a networkx graph object.

rasterize(resolution=50, bounding_square=None, warn=True, virtual_atoms=True, rotate=0, virtual_residues=True)[source]

Rasterize the projection to a square image of the given resolution. Uses the Bresenham algorithm for line rasterization.

Parameters:
  • resolution – The number of pixels in each direction.
  • bounding_square – Rasterize onto the given square. If None, automatically get a bounding_square that shows the whole projection
  • warn – If True, raise a warning if parts of the projection are not inside the given bounding square.
  • virtual_atoms – If True, virtual atoms are also rasterized.
  • rot – The in-plane rotation in degrees, applied before rotation.
Returns:

A tuple (np.array, float). The first value is a resolution x resolution numpy 2D array. The values are floats from 0.0 (black) to 1.0 (white). This array can be directly plotted using matplotlib: pyplot.imshow(array, cmap='gray', interpolation='none') The second value is the length of one pixle in angstrom.

rotate(angle)[source]

Rotate the projection in place around the origin (0,0)

Parameters:angle – Rotation angle in degrees.
vres_iterator
forgi.projection.projection2d.bresenham(start, end)[source]

Rasterize a line from start to end onto a grid with grid-width 1.

Parameters:
  • start – A sequence of length 2, containing the x,y coordinates of the start of the line
  • start – A sequence of length 2, containing the x,y coordinates of the end of the line
Returns:

A list of tuples (x,y), where x and y are integers.

forgi.projection.projection2d.crop_coordinates_to_bounds(a, num_cells)[source]
A:an array of x,y coordinates (modified in place)
forgi.projection.projection2d.diameter(Points)[source]

Given a list of 2d points, returns the pair that’s farthest apart.

forgi.projection.projection2d.hulls(Points)[source]

Graham scan to find upper and lower convex hulls of a set of 2d points.

forgi.projection.projection2d.log = <logging_exceptions.ExlogLogger object>

This module uses code by David Eppstein (http://code.activestate.com/recipes/117225-convex-hull-and-diameter-of-2d-point-sets/) under the PSF license and code by Syrtis Major (c)2014-2015 (http://stackoverflow.com/a/24567352/5069869) under the BSD 3-Clause license and potentially copyrighted code from matplotlib under the PSF license (http://matplotlib.org/examples/api/line_with_text.html).

forgi.projection.projection2d.orientation(p, q, r)[source]

Return positive if p-q-r are clockwise, neg if ccw, zero if colinear.

forgi.projection.projection2d.profile(x)[source]
forgi.projection.projection2d.rasterized_2d_coordinates(points, angstrom_per_cell=10, origin=array([0, 0]), rotate=0)[source]
forgi.projection.projection2d.rotate2D(vector, cosPhi, sinPhi)[source]
forgi.projection.projection2d.rotatingCalipers(Points)[source]

Given a list of 2d points, finds all ways of sandwiching the points between two parallel lines that touch one point each, and yields the sequence of pairs of points touched by each pair of lines.

forgi.projection.projection2d.to_grayscale(im)[source]

Convert an np array image from RGB to grayscale

forgi.projection.projection2d.to_rgb(im)[source]

Convert an np array image from grayscale to RGB

«  forgi.projection.hausdorff module   ::   Contents   ::   forgi.threedee package  »