stacker#
The Python Functions that control Command Line Options
This module contains the Python routines called when StACKER is run in the command line.
Functions
Disable printing to standard output |
|
Executes the residue movement routine to generate Bottaro values for a trajectory. |
|
|
Combines two 2D numpy arrays (frames_A and frames_B) into a single array. |
Compares pi-stacking events between two trajectories. |
|
Converts a parsed command to use to the correct subroutine and runs the routine |
|
Enable printing to standard output |
|
Executes the routine for filtering an input trajectory file and converting it to a PDB file. |
|
Calculates the distance between the centers of mass of two specified residues. |
|
Reads the user's passed in command line and runs the command |
|
Identifies residue pairs with the most π-stacking interactions. |
|
Runs the System Stacking Fingerprint (SSF) routine to generate a single SSF for a specified frame or range of frames. |
Exceptions
No Atomnames to subset to were provided |
|
No Frames present in trajectory |
|
Specified command line routine was invalid |
|
No Residues to subset to were provided |
- run_python_command()[source]#
Reads the user’s passed in command line and runs the command
Reads the command line input, runs the associated command with the added flags.
- convert_to_python_command()[source]#
Converts a parsed command to use to the correct subroutine and runs the routine
Converts the specified script to a python command and runs it with the associated inputs based on the flags.
- filter_traj_routine()[source]#
Executes the routine for filtering an input trajectory file and converting it to a PDB file.
This function uses provided command-line arguments to call filter_traj_to_pdb() with specified inputs, such as the trajectory file, topology file, output file location, residues to retain, and atom names to retain.
- Raises:
- ResEmpty
If the –residues argument is not provided.
- AtomEmpty
If the –atom_names argument is not provided.
See also
filter_traj_to_pdb
The core function that performs the trajectory filtering and PDB generation.
Notes
The args.residues argument should be a list of residue indices to retain.
The args.atom_names argument should be a comma-separated string of atom names to retain.
Output directories are created automatically if they do not exist.
Examples
Command-line usage with sample arguments:
$ stacker -s filter_traj -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -o testing/command_line_tests/filter/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb -r 426,427 -a C2,C4,C6
- bottaro_routine()[source]#
Executes the residue movement routine to generate Bottaro values for a trajectory.
This routine calculates the r, rho, and theta values for each frame of a PDB trajectory between two specified residues and stores them in a CSV file. Optionally, it visualizes the movement using heatmaps or scatter plots.
- Parameters:
- None
The routine relies on global args, which must be set via command-line arguments or equivalent argument parsing.
- Raises:
- AtomEmpty
If required atom names or residue indices are not provided.
See also
filter_traj_to_pdb
Filters trajectory and generates an intermediate PDB.
write_bottaro_to_csv
Writes residue movement calculations to a CSV.
visualize_two_residue_movement_heatmap
Creates a heatmap for residue movement.
visualize_two_residue_movement_scatterplot
Creates a scatter plot for residue movement.
Notes
args.trajectory and args.topology are required for generating an intermediate PDB if one is not already provided.
Supports heatmap and scatter plot visualizations for residue movement.
Removes intermediate files (pdb and csv) if –no_inter flag is specified.
Examples
Command-line usage:
$ stacker -s bottaro -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -pdb testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd.pdb -o testing/command_line_tests/bottaro/tUAG_aCUA_+1GCU_GC_plot.csv -p 426 -v 427 -pa C2,C4,C6 -va C2,C4,C6 -pt scatter $ stacker -s bottaro -pdb testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat_mdcrd_3200frames.pdb -o testing/command_line_tests/bottaro/tUAG_aCUA_+1GCU_GC_plot_3200frames.csv -p 426 -v 427 -pa C2,C4,C6 -va C2,C4,C6 -pt heat
- res_distance_routine()[source]#
Calculates the distance between the centers of mass of two specified residues.
This routine filters the trajectory based on the specified residues and atom names, and calculates the distances for either a single frame or multiple frames if bootstrapping is enabled.
- Parameters:
- None
The routine relies on global args for all input values.
- Raises:
- ResEmpty
If fewer or more than two residues are specified, or if no residues are provided.
- AtomEmpty
If no atom names are provided.
See also
filter_traj
Filters the trajectory for the specified residues and atom names.
calculate_residue_distance
Calculates the distance between residues for a given frame.
Notes
The args.residues argument must contain exactly two residues.
The args.atom_names argument is required to specify the atoms involved in the calculation.
If the –bootstrap argument is provided, distances are calculated across a specified number of randomly sampled frames.
If no bootstrap is performed, the distance is calculated for a single frame.
Examples
Command-line usage:
$ stacker -s res_distance -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -f 2 --residues 426,427 --atom_names C2,C4,C6
- system_routine()[source]#
Runs the System Stacking Fingerprint (SSF) routine to generate a single SSF for a specified frame or range of frames.
This routine processes a trajectory to create SSFs based on inter-residue distances and outputs either visualizations or raw data for further analysis. The user can provide input fingerprints or generate them from trajectory files.
- Parameters:
- None
The routine relies on global args for all input values.
- Raises:
- FileNotFoundError
If the specified trajectory or topology file does not exist.
- ValueError
If the provided input files have incompatible dimensions or unexpected formatting.
See also
filter_traj
Filters the trajectory for the specified residues.
get_residue_distance_for_trajectory
Computes inter-residue distances for a trajectory.
combine_frames
Combines two SSFs into a single array for comparison.
display_arrays_as_video
Visualizes SSF arrays and saves them as video files.
Notes
If args.input is provided, the routine will use pre-calculated fingerprint data.
If args.residues is specified, only those residues are considered in the calculation.
SSFs can be created for individual frames, a list of frames, or all frames in the trajectory.
If args.input_B is provided, two SSFs are combined for analysis.
Examples
Command-line usage:
$ stacker -s system -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -r 90-215 -fl 1-2 $ stacker -s ssf -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -r 90-215 -fl 1-2 -g 10 -o testing/command_line_tests/pairwise/5JUP_N2_tUAG_aCUA_+1GCU_nowat_pairwise_avg_1to2.png -d testing/command_line_tests/pairwise/5JUP_N2_tUAG_aCUA_+1GCU_data_1to2.txt.gz
- combine_frames(frames_A, frames_B)[source]#
Combines two 2D numpy arrays (frames_A and frames_B) into a single array.
This function takes two 2D numpy arrays of the same shape and combines them into a new array. The upper triangular part (excluding the diagonal) of the resulting array is filled with the corresponding elements from frames_A, while the lower triangular part (including the diagonal) is filled with the corresponding elements from frames_B.
- Parameters:
- frames_Anumpy.ndarray
A 2D numpy array.
- frames_Bnumpy.ndarray
A 2D numpy array of the same shape as frames_A.
- Returns:
- numpy.ndarray
A new 2D numpy array with combined elements from frames_A and frames_B.
Examples
>>> import stacker as st >>> import numpy as np >>> frames_A = np.array([[1, 2, 3], [4, 5, 6], [7, 8, 9]]) >>> frames_B = np.array([[9, 8, 7], [6, 5, 4], [3, 2, 1]]) >>> st.combine_frames(frames_A, frames_B) array([[9., 2., 3.], [6., 5., 6.],
[3., 2., 1.]])
- stack_events_routine()[source]#
Identifies residue pairs with the most π-stacking interactions.
This routine analyzes a molecular dynamics trajectory and identifies residue pairs with the smallest center of geometry (COG) distances, typically indicating strong π-stacking. It outputs the top n_events residue pairings based on user input.
- Parameters:
- None
This function relies on global args for input values.
- Raises:
- FileNotFoundError
If the trajectory or topology files are not found.
- ValueError
If input parameters are invalid or incompatible with the data.
See also
filter_traj
Filters trajectory data based on residue selection.
get_top_stacking
Extracts top stacking events from a trajectory.
Examples
Command-line usage:
$ stacker -s stack_events -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -r 90-215 -f 1 -n 5 $ stacker -s stack_events -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -r 90-100 -fl 1-10 -n 5 $ stacker -s stack_events -trj testing/first10_5JUP_N2_tUAG_aCUA_+1GCU_nowat.mdcrd -top testing/5JUP_N2_tUAG_aCUA_+1GCU_nowat.prmtop -r 90-215 -n 5 -i testing/command_line_tests/pairwise/5JUP_N2_tUAG_aCUA_+1GCU_data_1to2.txt
- compare_routine()[source]#
Compares pi-stacking events between two trajectories.
This routine analyzes two trajectory-derived stacking event files and identifies residue pairs with the largest change in average distances between the two systems. It processes the data, merges it, calculates discrepancies, and outputs the sorted results.
- Parameters:
- None
This function relies on global args for input values.
- Raises:
- FileNotFoundError
If one or both input files are not found.
- ValueError
If the input files have incompatible formats or lack the required headers.
See also
find_row_with_header
Identifies the header row in a data file.
preprocess_df
Prepares the input dataframes for comparison.
pd.merge
Merges the two dataframes for discrepancy calculations.
Notes
The input files must contain tab-separated values with the header Row Column Value.
Residue pairings are sorted alphabetically before comparison to ensure consistent results.
Examples
Command-line usage:
$ stacker -s compare -A /home66/esakkas/STACKER/SCRIPTS/slurmLogs_fingerprint/out_fingerprint_2418986 -B /home66/esakkas/STACKER/SCRIPTS/slurmLogs_fingerprint/out_fingerprint_2418997 -SA _tUAG_aCUA_+1GCU -SB _tUAG_aCUA_+1CGU
- exception InvalidRoutine[source]#
Bases:
Exception
Specified command line routine was invalid
This is raised when the stacker -s ROUTINE was given an invalid ROUTINE.
Examples
Command Line:
$ stacker -s blah
- __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.
- exception ResEmpty[source]#
Bases:
Exception
No Residues to subset to were provided
This is raised when a routine is given an invalid number of residues to subset to.
See also
fiter_traj_routine
bottaro_routine
Examples
Command Line:
$ stacker -s filter_traj -r -a C2,C4,C6
- __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.
- exception AtomEmpty[source]#
Bases:
Exception
No Atomnames to subset to were provided
This is raised when a routine is given an invalid number of atomnames to subset to.
See also
fiter_traj_routine
bottaro_routine
Examples
Command Line:
$ stacker -s filter_traj -a -r 426,427
- __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.
- exception FrameEmpty[source]#
Bases:
Exception
No Frames present in trajectory
- __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.
Modules
Filter Trajectory Files |
|
Run K Means Clustering and Principal Component Analysis |
|
Create and Analyze System Stacking Fingerprints (SSFs) |
|
Create and Analyze Pairwise Stacking Fingerprints (PSFs) |
|
Vector class used for SSF and PSF Calculation |
|
Visualize the SSFs, PSFs, and other analyses |