stacker.pairwise_distance module#
Create and Analyze System Stacking Fingerprints (SSFs)
This module contains the functions to create and analyze System Stacking Fingerprints. It allows the user to generate pairwise distance data across a trajectory and compare these matrices across different trajectories.
- calculate_residue_distance(trj, res1, res2, res1_atoms=('C2', 'C4', 'C6'), res2_atoms=('C2', 'C4', 'C6'), frame=1)[source]#
Calculates the vector between two residues with x, y, z units in Angstroms.
Calculates the distance between the center of two residues. The center is defined by the average x, y, z position of three passed atoms for each residue (typically every other carbon on the 6-carbon ring of the nucleotide base).
- Parameters:
- trjmd.Trajectory
Single frame trajectory.
- res1int
1-indexed residue number of the first residue (PDB Column 5).
- res2int
1-indexed residue number of the second residue (PDB Column 5).
- res1_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 1.
- res2_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 2.
- frameint, default=1
1-indexed frame number of trajectory to calculate the distance.
- Returns:
- distance_res12Vector
Vector from the center of geometry of residue 1 to residue 2.
See also
get_residue_distance_for_frame
Calculates pairwise distances between all residues in a given frame.
Examples
>>> import stacker as st >>> filtered_traj = st.filter_traj('testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... residues = {426,427}, ... atoms = {'C2','C4','C6'}) WARNING: Residue Indices are expected to be 1-indexed Reading trajectory... Reading topology... Filtering trajectory... WARNING: Output filtered traj atom, residue, and chain indices are zero-indexed >>> distance_vec = st.calculate_residue_distance( ... trajectory=filtered_traj, ... res1_num=426, ... res2_num=427, ... res1_atoms=("C2", "C4", "C6"), ... res2_atoms=("C2", "C4", "C6"), ... frame=1 ... ) >>> distance_vec.magnitude() 7.5253396
- get_residue_distance_for_frame(trj, frame, res1_atoms=('C2', 'C4', 'C6'), res2_atoms=('C2', 'C4', 'C6'), write_output=True)[source]#
Calculates System Stacking Fingerprint (SSF) between all residues in a given frame.
Calculates the distances between all pairs of residues in a given frame of a trajectory. Outputs as a square matrix with all residues on each side. This is the data behind a System Stacking Fingerprint (SSF)
- Parameters:
- trjmd.Trajectory
Trajectory to analyze (must have a topology).
- frameint
1-indexed frame to analyze.
- res1_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 1.
- res2_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 2.
- write_outputbool, default=True
If True, displays a loading screen to standard output.
- Returns:
- pairwise_distancesnp.typing.ArrayLike
Matrix where position (i, j) represents the distance from residue i to residue j.
See also
get_residue_distance_for_trajectory
Calculates System Stacking Fingerprints (SSFs) for all residues across all frames of a trajectory
filter_traj
Filters an input trajectory to only the specified atoms and residues
mdtraj.load
Load a trajectory+topology file
Examples
>>> import stacker as st >>> filtered_traj = st.filter_traj('stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... residues = '2-5,13-16,23-31,46-51,65-76,88-104,122-141,164-175,184-198,288-289,401-415,420-430', ... atoms = {'C2','C4','C6'}) >>> ssf = st.get_residue_distance_for_frame(filtered_traj, frame = 2, write_output = False) >>> ssf.shape (127, 127)
- get_residue_distance_for_trajectory(trj, frames={}, res1_atoms=('C2', 'C4', 'C6'), res2_atoms=('C2', 'C4', 'C6'), threads=1, write_output=True)[source]#
Calculates System Stacking Fingerprints (SSFs) for all residues across all frames of a trajectory.
Calculates the distances between all pairs of residues for all frames of a trajectory. Outputs as a square matrix with all residues on each side. This is the data behind a System Stacking Fingerprint (SSF)
- Parameters:
- trjmd.Trajectory
Trajectory to analyze (must have a topology). Output of st.filter_traj() or mdtraj.load()
- framesarray_like or str
list of frame indices to analyze (1-indexed). Accepts smart-indexed str representing a list of frames (e.g ‘1-5,6,39-48’) If empty, uses all frames.
- res1_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 1.
- res2_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 2.
- threadsint, default=1
Number of threads to use for parallel processing.
- write_outputbool, default=True
If True, displays a loading screen to standard output. Do not use if threads > 1.
- Returns:
- ssf_per_framearray_like
List where pairwise_distances[f] is the output of get_residue_distance_for_frame(trajectory, f, res1_atoms, res2_atoms).
See also
system_stacking_fingerprints
Alias for this function
get_residue_distance_for_frame
Calculates System Stacking Fingerprint (SSF) between all residues in a given frame.
filter_traj
Filters an input trajectory to only the specified atoms and residues
mdtraj.load
Load a trajectory+topology file
display_arrays_as_video
Displays this data as an SSF.
Examples
>>> import stacker as st >>> filtered_traj = st.filter_traj('stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... residues = '2-5,13-16,23-31,46-51,65-76,88-104,122-141,164-175,184-198,288-289,401-415,420-430', ... atoms = {'C2','C4','C6'}) >>> ssfs = st.get_residue_distance_for_trajectory(filtered_traj, frames = '1-3', write_output = False) >>> ssfs.shape (3, 127, 127)
- system_stacking_fingerprints(trj, frames={}, res1_atoms=('C2', 'C4', 'C6'), res2_atoms=('C2', 'C4', 'C6'), threads=1, write_output=True)[source]#
Alias for get_residue_distance_for_trajectory().
Calculates System Stacking Fingerprints (SSFs) for all residues across all frames of a trajectory.
Calculates the distances between all pairs of residues for all frames of a trajectory. Outputs as a square matrix with all residues on each side. This is the data behind a System Stacking Fingerprint (SSF)
- Parameters:
- trjmd.Trajectory
Trajectory to analyze (must have a topology). Output of st.filter_traj() or mdtraj.load()
- framesarray_like or str
list of frame indices to analyze (1-indexed). Accepts smart-indexed str representing a list of frames (e.g ‘1-5,6,39-48’) If empty, uses all frames.
- res1_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 1.
- res2_atomstuple, default=(“C2”, “C4”, “C6”)
Atom names whose positions are averaged to find the center of residue 2.
- threadsint, default=1
Number of threads to use for parallel processing.
- write_outputbool, default=True
If True, displays a loading screen to standard output. Do not use if threads > 1.
- Returns:
- ssf_per_framearray_like
List where pairwise_distances[f] is the output of get_residue_distance_for_frame(trajectory, f, res1_atoms, res2_atoms).
See also
system_stacking_fingerprints
Alias for this function
get_residue_distance_for_frame
Calculates System Stacking Fingerprint (SSF) between all residues in a given frame.
filter_traj
Filters an input trajectory to only the specified atoms and residues
mdtraj.load
Load a trajectory+topology file
display_arrays_as_video
Displays this data as an SSF.
Examples
>>> import stacker as st >>> filtered_traj = st.filter_traj('stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... residues = '2-5,13-16,23-31,46-51,65-76,88-104,122-141,164-175,184-198,288-289,401-415,420-430', ... atoms = {'C2','C4','C6'}) >>> ssfs = st.get_residue_distance_for_trajectory(filtered_traj, frames = '1-3', write_output = False) >>> ssfs.shape (3, 127, 127)
- get_frame_average(ssfs)[source]#
Calculates an average System Stacking Fingerprint (SSF) across multiple SSFs
Used to calculate an average SSF across multiple frames of a trajectory. Can average the result of get_residue_distance_for_trajectory
- Parameters:
- ssfsnumpy.typing.ArrayLike
List or array of 2D NumPy arrays representing a pairwise distance matrix of an MD structure. All 2D NumPy arrays must be of the same dimenstions. Output of
get_residue_distance_for_trajectory()
- Returns:
- avg_framenumpy.typing.ArrayLike
A single 2D NumPy array representing a pairwise distance matrix where each position i,j is the average distance from residue i to j across all matrices in frames.
See also
get_residue_distance_for_trajectory
Calculates System Stacking Fingerprints (SSFs) for all residues across all frames of a trajectory
system_stacking_fingerprints
Alias for get_residue_distance_for_trajectory
Examples
>>> import stacker as st >>> filtered_traj = st.filter_traj('stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... residues = '2-5,13-16,23-31,46-51,65-76,88-104,122-141,164-175,184-198,288-289,401-415,420-430', ... atoms = {'C2','C4','C6'}) >>> ssfs = st.get_residue_distance_for_trajectory(filtered_traj, frames = '1-3', write_output = False) >>> avg_ssf = st.get_frame_average(ssfs) >>> avg_ssf.shape [FILL]
- get_top_stacking(trj, matrix, csv='', n_events=5, include_adjacent=False)[source]#
Returns top stacking residue pairs for a given System Stacking Fingerprint (SSF)
Given a trajectory and a SSF made from get_residue_distance_for_frame() or get_frame_average() prints the residue pairings with the strongest stacking events (ie. the residue pairings with center of geometry distance closest to 3.5Å).
- Parameters:
- trjmd.Trajectory
trajectory used to get the stacking fingerprint
- matrixtyping.ArrayLike
Single-frame SSF created by
get_residue_distance_for_frame()
orget_frame_average()
- csvstr, default = ‘’,
output filename of the tab-separated txt file to write data to. If empty, data printed to standard output
- n_eventsint, default = 5
maximum number of stacking events to display, if -1 display all residue pairings
- include_adjacentbool, default = False
True if adjacent residues should be included in the printed output
See also
get_residue_distance_for_frame
Calculates System Stacking Fingerprint (SSF) between all residues in a given frame.
get_residue_distance_for_trajectory
Calculates System Stacking Fingerprints (SSFs) for all residues across all frames of a trajectory
system_stacking_fingerprints
Alias for get_residue_distance_for_trajectory
Examples
We can calculate the stacking events for a single frame:
>>> import stacker as st >>> filtered_traj = st.filter_traj('stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... atoms = {'C2','C4','C6'}) >>> ssf = st.get_residue_distance_for_frame(filtered_traj, frame = 2, write_output = False) >>> ssf.shape (252,252) >>> st.get_top_stacking(filtered_traj, ssf) Row Column Value 197 195 3.50 420 413 3.51 94 127 3.51 93 130 3.53 117 108 3.38
Or we can get most residue pairs that had the most stacking across many frames of a trajectory:
>>> import stacker as st >>> filtered_traj = st.filter_traj('stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd', ... 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop', ... atoms = {'C2','C4','C6'}) >>> ssfs = st.get_residue_distance_for_trajectory(filtered_traj, frames = '1-3', write_output = False) >>> avg_ssf = st.get_frame_average(ssfs) >>> avg_ssf.shape (252, 252) >>> st.get_top_stacking(filtered_traj, ssf) Row Column Value 130 93 3.56 108 117 3.44 195 197 3.61 127 94 3.65 47 167 3.66
- increment_residue(residue)[source]#
Increments residue ID by 1
Useful when converting from mdtraj 0-index residue naming to 1-indexed
- Parameters:
- residuestr
The residue id given by trajectory.topology.residue(i)
- Returns:
- incremented_idstr
The residue id with the sequence number increased by 1
Examples
>>> increment_residue('G43') 'G44'
- load_ssfs(file)[source]#
Loads a list of SSFs created by
stacker -s ssf -d OUTFILE
Loads a list of SSFs where each element is an SSF from a given frame. This is
OUTFILE
when runningstacker -s ssf -d OUTFILE
and quickly provides saved SSF data rather than recalculating SSFs for a given trajectory.- Parameters:
- filestr
outfile from
stacker -s ssf -d OUTFILE
. Must be.txt
or.txt.gz
.
- Returns:
- ssfsnumpy.typing.ArrayLike
List or array of 2D NumPy arrays representing a pairwise distance matrix of an MD structure. All 2D NumPy arrays must be of the same dimenstions. Output of
get_residue_distance_for_trajectory()
See also
system_stacking_fingerprints
calculates
ssfs
rather than loading it like this function.get_residue_distance_for_trajectory
Alias for
system_stacking_fingerprints()
get_frame_average
calculates average SSF for a trajectory.
matrix
parameter is the output of this
- exception MultiFrameTraj[source]#
Bases:
Exception
A multi-frame trajectory is passed to a one-frame function
Raised if a multi-frame trajectory is passed to a function that only works on one trajectory (eg. calculate_residue_distance_vector())
- __init__(*args, **kwargs)#
- add_note(object, /)#
Exception.add_note(note) – add a note to the exception
- args#
- with_traceback(object, /)#
Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.