visualization#

Visualize the SSFs, PSFs, and other analyses

This module includes the functions used to visualize plots, including SSFs and PSFs. The data inputs to these plot functions are provided by the other modules.

Functions

create_axis_labels(res_indices[, tick_distance])

Designates the axis labels to use in the SSF plot.

create_parent_directories(outfile_prefix)

Creates necessary parent directories to write an outfile given a prefix

display_arrays_as_video(ssfs, res_indices[, ...])

Displays SSF data to output or writes SSF as a PNG

display_psf_heatmap(csv[, plot_outfile, ...])

Alias for visualize_two_residue_movement_heatmap().

display_psf_scatter(csv[, plot_outfile, ...])

Alias for visualize_two_residue_movement_scatterplot().

display_ssfs(numpy_arrays, res_indices[, ...])

Alias for display_arrays_as_video().

set_polar_grid(**kwargs)

Set up axes for PSF.

visualize_two_residue_movement_heatmap(csv)

Creates heatmap of two-residue movement relative to each other.

visualize_two_residue_movement_scatterplot(csv)

Creates scatterplot of two-residue movement relative to each other.

Exceptions

NoResidues

Raised if user tries to make SSF with a trajectory of <1 residue

create_parent_directories(outfile_prefix)[source]#

Creates necessary parent directories to write an outfile given a prefix

Parameters:
outfile_prefixstr

The filepath of the output file, including the path where the file will be saved.

Examples

>>> create_parent_directories('/path/to/output/file.txt')
# This will create the directories '/path/to/output' if they do not exist.
>>> create_parent_directories('output/file.txt')
# This will create the directory 'output' if it does not exist.
create_axis_labels(res_indices, tick_distance=10)[source]#

Designates the axis labels to use in the SSF plot.

Helper function for visualizing SSFs. Returns the x-axis tick positions and labels to use on the ticks based on the residues used in a specific SSF analysis. Meant to be used when many disjoint sets of residue indices are used. Ticks will be present every tick_distance residues in a collection of adjacent residues, and a tick will exist at both ends of any consecutive residue sequence.

Parameters:
res_indiceslist

The list of residue indices used in the pairwise analysis. Parameter residue_desired passed to filter_traj()

tick_distanceint, default = 10

Distance between ticks in blocks of consecutive residues.

Returns:
tick_locationslist

List of tick positions (0-based) to place labels on in the axes.

tick_labelslist

List of labels to place at the adjacent tick locations.

See also

filter_traj

Filters an input trajectory to desired residues

Examples

Residues 0-12,98-100 were used. The SSF will label 0,10,12,98,100, provided in the second returned list. The first returned list gives the positions on the axes to place each label.

>>> create_axis_labels([0,1,2,3,4,5,6,7,8,9,10,11,12,98,99,100])
[0, 10, 12, 13, 15], [0, 10, 12, 98, 100]
>>> create_axis_labels([94,95,96,97,98,99,100,408,409,410,411,412,413,414,415,416,417,418,419,420,421,422,423,424,425,426,427,428])
[0, 6, 7, 17, 27], [94, 100, 408, 418, 428]
display_arrays_as_video(ssfs, res_indices, seconds_per_frame=10, tick_distance=10, outfile_prefix='', scale_limits=(0, 7), outfile='', scale_style='bellcurve', xy_line=True, **kwargs)[source]#

Displays SSF data to output or writes SSF as a PNG

Visualizes the data for an SSF for a trajectory or a single frame. Takes an SSF array outputted from get_residue_distance_for_frame, get_residue_distance_for_trajectory, or system_stacking_fingerprints and treats them as frames in a video, filling in a grid at position i, j by the value at i, j in the array.

Parameters:
ssfsarray_like

List or array of 2D NumPy arrays representing SSFs, output of system_stacking_fingerprints

res_indiceslist or str

The list of residue indices used in the pairwise analysis. Parameter residue_desired passed to filter_traj() Accepts smart-indexed str representing a list of residues (e.g ‘1-5,6,39-48’)

seconds_per_frameint, default = 10

Number of seconds to display each matrix for.

tick_distanceint, default = 10

Distance between ticks in blocks of consecutive residues.

outfilestr

Image output filepath to write a single SSF to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

outfile_prefixstr

Prefix for image filepath to write multiple frames to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

scale_limitstuple, default = (0, 7)

Limits of the color scale.

scale_style{‘bellcurve’, ‘gradient’}, default = ‘bellcurve’

Style of color scale.

xy_linebool, default = True

Draw x = y line to separate matrix halves.

**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 8

    Width of the figure in inches.

  • fig_heightfloat, default = 8

    Height of the figure in inches.

  • title_fontsizeint, default = fontsize

    Font size for the title.

  • legend_fontsizeint, default = fontsize

    Font size for the legend.

  • cb_fontsizeint, default = fontsize

    Font size for the colorbar tick labels.

  • xaxis_fontsizeint, default = fontsize

    Font size for the x-axis tick labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the y-axis tick labels.

Returns:
None

Displays video of NumPy arrays.

See also

create_axis_labels

Designates the axis labels to use in the SSF plot.

get_residue_distance_for_frame

Calculates System Stacking Fingerprint (SSF) between all residues in a given frame.

get_residue_distance_for_trajectory

get SSF data for all frames of a trajectory

system_stacking_fingerprints

Alias for this get_residue_distance_for_trajectory

display_ssfs

Alias for this function.

Examples

>>> import stacker as st
>>> import mdtraj as md
>>> trajectory_file = 'stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd'
>>> topology_file = 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop'
>>> trj = md.load(trajectory_file, top = topology_file)
>>> residue_selection_query = 'resi 90 to 215'
>>> frames_to_include = [1,2,3,4,5]
>>> trj_sub = trj.atom_slice(trj.top.select(residue_selection_query))
>>> resSeqs = [res.resSeq for res in trj_sub.topology.residues]
>>> frames = st.get_residue_distance_for_trajectory(trj_sub, frames_to_include, threads = 5)
Frame 2 done.
Frame 3 done.
Frame 1 done.
Frame 5 done.
Frame 4 done.
>>> st.display_arrays_as_video([st.get_frame_average(frames)], resSeqs, seconds_per_frame=10)
# Displays SSF for each frame of this trajectory to standard output
display_ssfs(numpy_arrays, res_indices, seconds_per_frame=10, tick_distance=10, outfile_prefix='', scale_limits=(0, 7), outfile='', scale_style='bellcurve', xy_line=True, **kwargs)[source]#

Alias for display_arrays_as_video().

display_arrays_as_video(

ssfs, res_indices, seconds_per_frame = 10, tick_distance = 10, outfile_prefix = ‘’, scale_limits = (0,7), outfile = ‘’, scale_style = ‘bellcurve’, xy_line = True, **kwargs

)

Displays SSF data to output or writes SSF as a PNG

Visualizes the data for an SSF for a trajectory or a single frame. Takes an SSF array outputted from get_residue_distance_for_frame, get_residue_distance_for_trajectory, or system_stacking_fingerprints and treats them as frames in a video, filling in a grid at position i, j by the value at i, j in the array.

Parameters:
ssfsarray_like

List or array of 2D NumPy arrays representing SSFs, output of system_stacking_fingerprints

res_indiceslist or str

The list of residue indices used in the pairwise analysis. Parameter residue_desired passed to filter_traj() Accepts smart-indexed str representing a list of residues (e.g ‘1-5,6,39-48’)

seconds_per_frameint, default = 10

Number of seconds to display each matrix for.

tick_distanceint, default = 10

Distance between ticks in blocks of consecutive residues.

outfilestr

Image output filepath to write a single SSF to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

outfile_prefixstr

Prefix for image filepath to write multiple frames to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

scale_limitstuple, default = (0, 7)

Limits of the color scale.

scale_style{‘bellcurve’, ‘gradient’}, default = ‘bellcurve’

Style of color scale.

xy_linebool, default = True

Draw x = y line to separate matrix halves.

**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 8

    Width of the figure in inches.

  • fig_heightfloat, default = 8

    Height of the figure in inches.

  • title_fontsizeint, default = fontsize

    Font size for the title.

  • legend_fontsizeint, default = fontsize

    Font size for the legend.

  • cb_fontsizeint, default = fontsize

    Font size for the colorbar tick labels.

  • xaxis_fontsizeint, default = fontsize

    Font size for the x-axis tick labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the y-axis tick labels.

Returns:
None

Displays video of NumPy arrays.

See also

create_axis_labels

Designates the axis labels to use in the SSF plot.

get_residue_distance_for_frame

Calculates System Stacking Fingerprint (SSF) between all residues in a given frame.

get_residue_distance_for_trajectory

get SSF data for all frames of a trajectory

system_stacking_fingerprints

Alias for this get_residue_distance_for_trajectory

display_ssfs

Alias for this function.

Examples

>>> import stacker as st
>>> import mdtraj as md
>>> trajectory_file = 'stacker/testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd'
>>> topology_file = 'stacker/testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop'
>>> trj = md.load(trajectory_file, top = topology_file)
>>> residue_selection_query = 'resi 90 to 215'
>>> frames_to_include = [1,2,3,4,5]
>>> trj_sub = trj.atom_slice(trj.top.select(residue_selection_query))
>>> resSeqs = [res.resSeq for res in trj_sub.topology.residues]
>>> frames = st.get_residue_distance_for_trajectory(trj_sub, frames_to_include, threads = 5)
Frame 2 done.
Frame 3 done.
Frame 1 done.
Frame 5 done.
Frame 4 done.
>>> st.display_arrays_as_video([st.get_frame_average(frames)], resSeqs, seconds_per_frame=10)
# Displays SSF for each frame of this trajectory to standard output
set_polar_grid(**kwargs)[source]#

Set up axes for PSF.

Creates polar plot background for two-residue movement comparison with theta 0 to 360, a radial maximum of 15 Angstroms, and a visualization of the perspective residue at the center.

Parameters:
**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 5

    Width of the figure in inches.

  • fig_heightfloat, default = 5

    Height of the figure in inches.

  • title_fontsizeint, default = fontsize

    Font size for the title.

  • legend_fontsizeint, default = fontsize

    Font size for the legend.

  • cb_fontsizeint, default = fontsize

    Font size for the colorbar labels.

  • xaxis_fontsizeint, default = fontsize

    Font size for the x-axis labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the y-axis labels.

Returns:
axmatplotlib.projections.polar.PolarAxes

Axis object for the created polar plot.

visualize_two_residue_movement_scatterplot(csv, plot_outfile='', frame_list={}, **kwargs)[source]#

Creates scatterplot of two-residue movement relative to each other.

Takes the data created in residue_movement and visualizes it as a polar coordinate scatterplot similar to the Figure D link in Proposal Feature 4.

Parameters:
csvstr

Filepath to csv file containing data on the movement of two residues relative to each other (r, rho, and theta values). Created in residue_movement.

plot_outfilestr

Filepath of the image file to write to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

frame_listset, default = {}

Set of frames to use in csv, if empty use all frames.

**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 5

    Width of the figure in inches.

  • fig_heightfloat, default = 5

    Height of the figure in inches.

  • xaxis_fontsizeint, default = fontsize

    Font size for the theta-axis labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the rho-axis labels.

Returns:
None

See also

write_bottaro_to_csv

Creates CSV file that is inputted here

Examples

>>> import stacker as st
>>> trajectory_file = 'testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd'
>>> topology_file = 'testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop'
>>> pdb_filename = 'testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb'
>>> output_csv_name = "testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv"
>>> perspective_residue = 426 # 1-indexed
>>> viewed_residue = 427 # 1-indexed
>>> st.filter_traj_to_pdb(trj_file=trajectory_file, top_file=topology_file, pdb=pdb_filename,
...                        residues={perspective_residue,viewed_residue}, 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
WARNING: Output file atom, residue, and chain indices are zero-indexed
Filtered trajectory written to:  testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb
>>> st.write_bottaro_to_csv(pdb_filename, output_csv_name, perspective_residue_num=perspective_residue, viewed_residue_num=viewed_residue)
Output values written to testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv
>>> st.visualize_two_residue_movement_scatterplot('testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv', 
...                                                 plot_outfile='testing/script_tests/visualization/tUAG_aCUA_+1GCU_GC_plot_10frames_scatter.png')
visualize_two_residue_movement_heatmap(csv, plot_outfile='', frame_list={}, **kwargs)[source]#

Creates heatmap of two-residue movement relative to each other.

2D shaded contour plot of the density of points in the visualize_two_residue_movement_scatterplot() scatterplot.

Parameters:
csvstr

Filepath to csv file containing data on the movement of two residues relative to each other (r, rho, and theta values). Created in residue_movement.

plot_outfilestr

Filepath of the image file to write to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

frame_listset, default = {}

Set of frames to use in csv, if empty use all frames.

**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 5

    Width of the figure in inches.

  • fig_heightfloat, default = 5

    Height of the figure in inches.

  • xaxis_fontsizeint, default = fontsize

    Font size for the theta-axis labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the rho-axis labels.

  • cb_fontsizeint, default = 10

    Font size for the colorbar tick labels.

Returns:
None

See also

write_bottaro_to_csv

Creates CSV file that is inputted here

Examples

>>> import stacker as st
>>> trajectory_file = 'testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd'
>>> topology_file = 'testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop'
>>> pdb_filename = 'testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb'
>>> output_csv_name = "testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv"
>>> perspective_residue = 426 # 1-indexed
>>> viewed_residue = 427 # 1-indexed
>>> st.filter_traj_to_pdb(trj_file=trajectory_file, top_file=topology_file, pdb=pdb_filename,
...                        residues={perspective_residue,viewed_residue}, 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
WARNING: Output file atom, residue, and chain indices are zero-indexed
Filtered trajectory written to:  testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb
>>> st.write_bottaro_to_csv(pdb_filename, output_csv_name, pers_res=perspective_residue, view_res=viewed_residue)
Output values written to testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv
>>> st.visualize_two_residue_movement_heatmap('testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv', 
...                                                 plot_outfile='testing/script_tests/visualization/tUAG_aCUA_+1GCU_GC_plot_10frames_heat.png')
display_psf_scatter(csv, plot_outfile='', frame_list={}, **kwargs)[source]#

Alias for visualize_two_residue_movement_scatterplot().

Creates scatterplot of two-residue movement relative to each other.

Takes the data created in residue_movement and visualizes it as a polar coordinate scatterplot similar to the Figure D link in Proposal Feature 4.

Parameters:
csvstr

Filepath to csv file containing data on the movement of two residues relative to each other (r, rho, and theta values). Created in residue_movement.

plot_outfilestr

Filepath of the image file to write to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

frame_listset, default = {}

Set of frames to use in csv, if empty use all frames.

**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 5

    Width of the figure in inches.

  • fig_heightfloat, default = 5

    Height of the figure in inches.

  • xaxis_fontsizeint, default = fontsize

    Font size for the theta-axis labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the rho-axis labels.

Returns:
None

See also

write_bottaro_to_csv

Creates CSV file that is inputted here

Examples

>>> import stacker as st
>>> trajectory_file = 'testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd'
>>> topology_file = 'testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop'
>>> pdb_filename = 'testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb'
>>> output_csv_name = "testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv"
>>> perspective_residue = 426 # 1-indexed
>>> viewed_residue = 427 # 1-indexed
>>> st.filter_traj_to_pdb(trj_file=trajectory_file, top_file=topology_file, pdb=pdb_filename,
...                        residues={perspective_residue,viewed_residue}, 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
WARNING: Output file atom, residue, and chain indices are zero-indexed
Filtered trajectory written to:  testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb
>>> st.write_bottaro_to_csv(pdb_filename, output_csv_name, perspective_residue_num=perspective_residue, viewed_residue_num=viewed_residue)
Output values written to testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv
>>> st.visualize_two_residue_movement_scatterplot('testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv', 
...                                                 plot_outfile='testing/script_tests/visualization/tUAG_aCUA_+1GCU_GC_plot_10frames_scatter.png')
display_psf_heatmap(csv, plot_outfile='', frame_list={}, **kwargs)[source]#

Alias for visualize_two_residue_movement_heatmap().

Creates heatmap of two-residue movement relative to each other.

2D shaded contour plot of the density of points in the visualize_two_residue_movement_scatterplot() scatterplot.

Parameters:
csvstr

Filepath to csv file containing data on the movement of two residues relative to each other (r, rho, and theta values). Created in residue_movement.

plot_outfilestr

Filepath of the image file to write to. Format inferred from file extension. png, pdf, ps, eps, and svg supported.

frame_listset, default = {}

Set of frames to use in csv, if empty use all frames.

**kwargsdict, optional

Additional keyword arguments to customize the plot:

  • fontsizeint, default = 10

    Font size for all text elements.

  • fig_widthfloat, default = 5

    Width of the figure in inches.

  • fig_heightfloat, default = 5

    Height of the figure in inches.

  • xaxis_fontsizeint, default = fontsize

    Font size for the theta-axis labels.

  • yaxis_fontsizeint, default = fontsize

    Font size for the rho-axis labels.

  • cb_fontsizeint, default = 10

    Font size for the colorbar tick labels.

Returns:
None

See also

write_bottaro_to_csv

Creates CSV file that is inputted here

Examples

>>> import stacker as st
>>> trajectory_file = 'testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd'
>>> topology_file = 'testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop'
>>> pdb_filename = 'testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb'
>>> output_csv_name = "testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv"
>>> perspective_residue = 426 # 1-indexed
>>> viewed_residue = 427 # 1-indexed
>>> st.filter_traj_to_pdb(trj_file=trajectory_file, top_file=topology_file, pdb=pdb_filename,
...                        residues={perspective_residue,viewed_residue}, 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
WARNING: Output file atom, residue, and chain indices are zero-indexed
Filtered trajectory written to:  testing/script_tests/residue_movement/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb
>>> st.write_bottaro_to_csv(pdb_filename, output_csv_name, pers_res=perspective_residue, view_res=viewed_residue)
Output values written to testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv
>>> st.visualize_two_residue_movement_heatmap('testing/script_tests/residue_movement/tUAG_aCUA_+1GCU_GC_plot.csv', 
...                                                 plot_outfile='testing/script_tests/visualization/tUAG_aCUA_+1GCU_GC_plot_10frames_heat.png')
exception NoResidues[source]#

Bases: Exception

Raised if user tries to make SSF with a trajectory of <1 residue

__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.