Try pytraj
online:
Submodules
pytraj.
distance
(traj=None, mask='', frame_indices=None, dtype='ndarray', top=None, image=False, n_frames=None)¶compute distance between two maskes
Parameters: | traj : Trajectory-like, list of Trajectory, list of Frames mask : str or a list of string or a 2D array-like of integers frame_indices : array-like, optional, default None dtype : return type, default ‘ndarray’ top : Topology, optional image : bool, default False n_frames : int, optional, default None
|
---|---|
Returns: | 1D ndarray if mask is a string 2D ndarray, shape (n_atom_pairs, n_frames) if mask is a list of strings or an array |
Notes
Be careful with Topology. If your topology has Box info but your traj does not, you would get weird output ([0.0, ...]). Make sure to use image=False in this method or set_nobox for Topology.
Examples
>>> import pytraj as pt
>>> # calculate distance for two atoms, using amber mask
>>> traj = pt.datafiles.load_tz2_ortho()
>>> dist = pt.distance(traj, '@1 @3')
>>> # calculate distance for two groups of atoms, using amber mask
>>> dist = pt.distance(traj, '@1,37,8 @2,4,6')
>>> # calculate distance between two residues, using amber mask
>>> dist = pt.distance(traj, ':1 :10')
>>> # calculate multiple distances between two residues, using amber mask
>>> # distance between residue 1 and 10, distance between residue 3 and 20
>>> # (when using atom string mask, index starts from 1)
>>> dist = pt.distance(traj, [':1 :10', ':3 :20'])
>>> # calculate distance for a series of atoms, using array for atom mask
>>> # distance between atom 1 and 5, distance between atom 4 and 10 (index starts from 0)
>>> dist = pt.distance(traj, [[1, 5], [4, 10]])
pytraj.
angle
(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None, *args, **kwargs)¶compute angle between two maskes
Parameters: | traj : Trajectory-like, list of Trajectory, list of Frames mask : str or array top : Topology, optional dtype : return type, defaul ‘ndarray’ |
---|---|
Returns: | 1D ndarray if mask is a string 2D ndarray, shape (n_atom_pairs, n_frames) if mask is a list of strings or an array |
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # calculate angle for three atoms, using amber mask
>>> pt.angle(traj, '@1 @3 @10')
array([ 98.06193365, 95.75979717, 105.59626997, 107.64079091,
94.93516228, 102.06028369, 109.3479057 , 110.11532478,
101.86104127, 110.48992512])
>>> # calculate angle for three groups of atoms, using amber mask
>>> angles = pt.angle(traj, '@1,37,8 @2,4,6 @5,20')
>>> # calculate angle between three residues, using amber mask
>>> angles = pt.angle(traj, ':1 :10 :20')
>>> # calculate multiple angles between three residues, using amber mask
>>> # angle between residue 1, 10, 20, angle between residue 3, 20, 30
>>> # (when using atom string mask, index starts from 1)
>>> angles = pt.angle(traj, [':1 :10 :20', ':3 :20 :30'])
>>> # calculate angle for a series of atoms, using array for atom mask
>>> # angle between atom 1, 5, 8, distance between atom 4, 10, 20 (index starts from 0)
>>> angles = pt.angle(traj, [[1, 5, 8], [4, 10, 20]])
pytraj.
rmsd
(traj=None, mask='', ref=0, ref_mask='', nofit=False, mass=False, update_coordinate=True, frame_indices=None, top=None, dtype='ndarray')¶compute rmsd
Parameters: | traj : Trajectory-like mask : str or 1D array-like of string or 1D or 2D array-like
ref : {Frame, int}, default=0 (first frame)
ref_mask: str, optional
nofit : bool, default False
mass : bool, default False
update_coordinate : bool, default True
frame_indices : int 1D array-like, default None
top : {Topology, str}, default None, optional dtype : return data type, default=’ndarray’ |
---|
Notes
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_trpcage()
>>> # all atoms, do fitting, using ref=traj[-3]
>>> data = pt.rmsd(traj, ref=-3)
>>> # rmsd for 3 maskes, do fitting, using ref=traj[0] (defaul)
>>> data = pt.rmsd(traj, mask=['@CA', '@C', ':3-18@CA'], dtype='dataset')
>>> # rmsd to first frame, use mass ':3-13' but do not perorm fitting
>>> data= pt.rmsd(traj, ref=traj[0], mask=':3-13', nofit=True)
>>> # use atom indices for mask
>>> data= pt.rmsd(traj, ref=traj[0], mask=range(40), nofit=True)
>>> # compute rmsd (and align) with reference having different atoms
>>> trpcage_traj = pt.datafiles.load_trpcage()[:]
>>> tz2_traj = pt.datafiles.load_tz2()[:1]
>>> data = pt.rmsd(trpcage_traj, mask='@1-10', ref=tz2_traj, ref_mask='@11-20')
>>> data
array([ 2.16203842, 2.28859396, 2.15817654, ..., 2.20767189,
2.30087764, 1.92654945])
Notes
-----
if ``traj`` is mutable and update_coordinate=True, its coordinates will be updated.
pytraj.
rmsd_perres
(traj=None, mask='', ref=0, mass=False, resrange=None, perres_mask=None, perres_center=False, perres_invert=False, frame_indices=None, top=None, dtype='dataset', **kwd)¶superpose traj
to ref
with mask, then calculate nofit rms for residues
in resrange with given perresmask
Returns: | out : pytraj.DatasetList, shape=(1+n_residues, n_frames)
|
---|
pytraj.
pairwise_rmsd
(traj=None, mask='', metric='rms', top=None, dtype='ndarray', mat_type='full', frame_indices=None)¶Calculate pairwise rmsd with different metrics.
Parameters: | traj : Trajectory-like or iterable object mask : mask
metric : {‘rms’, ‘dme’, ‘srmsd’, ‘nofit’}
top : Topology, optional, default=None dtype: ndarray
mat_type : str, {‘full’, ‘half’}
|
---|
Notes
Install libcpptraj
with openmp
to get benefit from parallel
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> arr = pt.pairwise_rmsd(traj(0, 1000, mask='@CA'))
>>> # calculate pairwise rmsd for all frames using CA atoms, use `dme` (distance RMSD)
>>> # convert to numpy array
>>> arr_np = pt.pairwise_rmsd(traj, "@CA", metric="dme", dtype='ndarray')
>>> # calculate pairwise rmsd for all frames using CA atoms, nofit for RMSD
>>> # convert to numpy array
>>> arr_np = pt.pairwise_rmsd(traj, "@CA", metric="nofit", dtype='ndarray')
>>> # calculate pairwise rmsd for all frames using CA atoms
>>> # use symmetry-corrected RMSD, convert to numpy array
>>> arr_np = pt.pairwise_rmsd(traj, "@CA", metric="srmsd", dtype='ndarray')
>>> # use different dtype
>>> arr_np = pt.pairwise_rmsd(traj, "@CA", metric="srmsd", dtype='dataset')
pytraj.
dssp
(traj=None, mask='', frame_indices=None, dtype='ndarray', simplified=False, top=None)¶return dssp profile for frame/traj
Parameters: | traj : Trajectory-like mask: str
frame_indices : {None, array-like}, default None, optional
dtype : str, default ‘ndarray’
simplified : bool, default False
|
---|---|
Returns: | out_0: ndarray, shape=(n_residues,)
out_1: ndarray, shape=(n_frames, n_residues)
out_2 : pytraj.DatasetList
|
Notes
Character | Integer | DSSP_Char | Secondary structure type |
---|---|---|---|
0 | 0 | ‘0’ | None |
b | 1 | ‘E’ | Parallel Beta-sheet |
B | 2 | ‘B’ | Anti-parallel Beta-sheet |
G | 3 | ‘G’ | 3-10 helix |
H | 4 | ‘H’ | Alpha helix |
I | 5 | ‘I’ | Pi (3-14) helix |
T | 6 | ‘T’ | Turn |
S | 7 | ‘S’ | Bend |
Simplified codes:
- 'H': include 'H', 'G', 'I' (helix)
- 'E': include 'E', 'B' (strand)
- 'C': include 'T', 'S' or '0' (coil)
Simplified codes will be mostly used for visualization in other packages.
Examples
>>> import pytraj as pt
>>> traj = pt.load_pdb_rcsb('1l2y')
>>> residues, ss, _ = pt.dssp(traj, ":2-10")
>>> residues
array(['LEU:2', 'TYR:3', 'ILE:4', 'GLN:5', 'TRP:6', 'LEU:7', 'LYS:8',
'ASP:9', 'GLY:10'],
dtype='<U6')
>>> ss
array([['0', 'H', 'H', ..., 'H', 'T', '0'],
['0', 'H', 'H', ..., 'H', 'T', '0'],
['0', 'H', 'H', ..., 'H', 'T', '0'],
...,
['0', 'H', 'H', ..., 'H', 'T', '0'],
['0', 'H', 'H', ..., 'H', 'H', '0'],
['0', 'H', 'H', ..., 'H', 'T', '0']],
dtype='<U1')
>>> residues, ss, _ = pt.dssp(traj, mask=range(100))
>>> traj = pt.fetch_pdb('1l2y')
>>> residues, ss, _ = pt.dssp(traj, simplified=True)
>>> ss[0].tolist() # first frame
['C', 'H', 'H', 'H', 'H', 'H', 'H', 'H', 'C', 'C', 'H', 'H', 'H', 'H', 'C', 'C', 'C', 'C', 'C', 'C']
pytraj.
multidihedral
(traj=None, dihedral_types=None, resrange=None, define_new_type=None, range360=False, dtype='dataset', top=None, frame_indices=None)¶perform dihedral search
Parameters: | traj : Trajectory-like object dihedral_types : dihedral type, default None
resrange : str | array-like
define_new_type : str
range360 : bool, default False
top : Topology | str, optional
|
---|---|
Returns: | pytraj.DatasetList (use values attribute to get raw numpy array) |
Notes
Dataset lables show residue number in 1-based index
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.multidihedral(traj)
>>> data = pt.multidihedral(traj, 'phi psi')
>>> data = pt.multidihedral(traj, resrange=range(8))
>>> data = pt.multidihedral(traj, range360=True)
>>> data = pt.multidihedral(traj, resrange='1,3,5')
>>> data = pt.multidihedral(traj, dihedral_types='phi psi')
>>> data = pt.multidihedral(traj, dihedral_types='phi psi', resrange='3-7')
>>> data = pt.multidihedral(traj, dihedral_types='phi psi', resrange=[3, 4, 8])
pytraj.
bfactors
(traj=None, mask='', byres=True, top=None, dtype='ndarray', frame_indices=None)¶calculate pseudo bfactor
Parameters: | traj: Trajectory-like mask: str, mask |
---|---|
Returns: | if dtype is ‘ndarray’ (default), return a numpy array with shape=(n_atoms/n_residues, 2) ([atom_or_residue_idx, value]) |
Notes
This is NOT getting bfactor from xray, but computing bfactor from simulation.
Examples
>>> import pytraj as pt
>>> from pytraj.testing import get_fn
>>> fn, tn = get_fn('tz2')
>>> traj = pt.load(fn, tn, mask='!:WAT')
>>> traj = pt.superpose(traj)
>>> bfactor = pt.bfactors(traj, byres=True)
pytraj.
radgyr
(traj=None, mask='', top=None, nomax=True, frame_indices=None, dtype='ndarray')¶compute radius of gyration
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.radgyr(traj, '@CA')
>>> data = pt.radgyr(traj, '!:WAT', nomax=False)
>>> data = pt.radgyr(traj, '@CA', frame_indices=[2, 4, 6])
pytraj.
molsurf
(traj=None, mask='', probe=1.4, offset=0.0, dtype='ndarray', frame_indices=None, top=None)¶calc molsurf
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> pt.molsurf(traj, '@CA')
array([ 458.51409637, 459.64784573, 456.54690793, 467.72939574,
462.45908781, 458.70327554, 454.40514806, 455.15015576,
468.70566447, 456.0058624 ])
>>> pt.molsurf(traj, '!:WAT')
array([ 1079.1395679 , 1090.79759341, 1069.65127413, 1096.0810919 ,
1091.65862234, 1091.68906298, 1085.53105392, 1069.22510187,
1079.70803583, 1075.8151414 ])
pytraj.
center_of_mass
(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None)¶compute center of mass
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> # compute center of mass residue 3 for first 2 frames.
array([[-0.661702 , 6.69124347, 3.35159413],
[ 0.5620708 , 7.82263042, -0.72707798]])
pytraj.
center_of_geometry
(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None)¶pytraj.
search_hbonds
(traj, mask='', solvent_donor=None, solvent_acceptor=None, distance=3.0, angle=135.0, image=False, series=True, options='', dtype='hbond', frame_indices=None, top=None)¶(combined with cpptraj doc) Searching for Hbond donors/acceptors in region specified by mask
.
Hydrogen bond is defined as A-HD, where A is acceptor heavy atom, H is hydrogen, D is
donor heavy atom. Hydrogen bond is formed when A to D distance < distance cutoff and A-H-D angle
> angle cutoff; if angle < 0 it is ignored.
Parameters: | traj : Trajectory-like mask : {str, 1D array-like}
solvent_donor : {None, str}, default None solvent_acceptor: {None, str}, deafult None
distance : float, default 3.0 (angstrom)
angle : float, 135.0 degree
dtype : return output’s type, default ‘hbond’ image : bool, default False series : bool, default True
options : str
|
---|---|
Returns: | out : DatasetHBond if series is True else return ‘DatasetList’ |
See also
to_amber_mask
Notes
Examples
>>> import pytraj as pt
>>> traj = pt.load_sample_data('tz2')
>>> # search hbond without including solvent
>>> data = pt.search_hbonds(traj, ':5,8')
>>> data
<pytraj.hbonds.DatasetHBond
donor_acceptor pairs : 2>
>>> data.donor_acceptor
['LYS8_O-GLU5_N-H', 'GLU5_O-LYS8_N-H']
>>> # get raw data, ndarray with shape=(n_hbonds+1, n_frames)
>>> # first array shows the total solute hbonds and other arrays shows
>>> # if hbond exists (1) or non-exists (0) for each frame
>>> data.values
array([[2, 2, 0, ..., 1, 1, 1],
[1, 1, 0, ..., 1, 1, 1],
[1, 1, 0, ..., 0, 0, 0]], dtype=int32)
>>> # search hbond including solvent
>>> hbonds = pt.search_hbonds(traj, ':5,8', solvent_donor=':WAT@O', solvent_acceptor=':WAT')
>>> hbonds
<pytraj.hbonds.DatasetHBond
donor_acceptor pairs : 8>
>>> hbonds.donor_acceptor
['LYS8_O-GLU5_N-H', 'GLU5_O-LYS8_N-H', 'LYS8_HZ1-V', 'LYS8_HZ2-V', 'GLU5_OE2-V', 'GLU5_O-V', 'GLU5_OE1-V', 'LYS8_HZ3-V']
>>> # 'GLU5_O-V' mean non-specific hbond between GLU5_O and solvent (:WAT in this case)
pytraj.
jcoupling
(traj=None, mask='', top=None, kfile=None, dtype='dataset', frame_indices=None)¶compute j-coupling
Parameters: | traj : any things that make frame_iter_master returning Frame object command : str, default “”
kfile : str, default None, optional
dtype : str, {‘dataset’, ...}, default ‘dataset’ |
---|
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> data = pt.calc_jcoupling(traj, ':1-12', kfile='data/Karplus.txt')
pytraj.
density
(traj, mask='*', density_type='number', delta=0.25, direction='z', dtype='dict')¶Compute density (number, mass, charge, electron) along a coordinate
Parameters: | traj : Trajectory-like mask : str or list of str, default ‘*’
density_type : str, {‘number’, ‘mass’, ‘charge’, ‘electron’}, default ‘number’ delta : float, default 0.25
direction : str, default ‘z’ dtype : str, default ‘dict’
|
---|---|
Returns: | out : dict of average density and std for each frame |
Notes
Syntax might be changed
Examples
>>> def func():
... import pytraj as pt
... fn = "data/DOPC.rst7"
... tn = "data/DOPC.parm7"
... traj = pt.load("data/DOPC.rst7", "data/DOPC.parm7")
... delta = ‘0.25’ ... density_type = ‘charge’ ... masks = [”:PC@P31”, ”:PC@N31”, ”:PC@C2”, ”:PC | :OL | :OL2”] ... density_dict = pt.density(traj, mask=masks, density_type=density_type, delta=delta) ... return density_dict >>> density_dict = func() # doctest: +SKIP
pytraj.
volume
(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None)¶compute volume
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> vol = pt.volume(traj, '@CA')
pytraj.
mindist
(traj=None, command='', top=None, dtype='ndarray', frame_indices=None)¶compute mindist
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> data = pt.mindist(traj, '@CA @H')
pytraj.
acorr
(data, dtype='ndarray', option='')¶compute autocorrelation
Parameters: | data : 1d array-like dtype: return type, default ‘ndarray’ covar : bool, default True option : str
|
---|
Notes
Same as autocorr in cpptraj
pytraj.
xcorr
(data0, data1, dtype='ndarray')¶compute cross correlation between two datasets
Parameters: | data0 and data1: 1D-array like dtype : return datatype, default ‘ndarray’ |
---|
Notes
Same as corr in cpptraj
pytraj.
timecorr
(vec0, vec1, order=2, tstep=1.0, tcorr=10000.0, norm=False, dtype='ndarray')¶compute time correlation.
Parameters: | vec0 : 2D array-like, shape=(n_frames, 3) vec1 : 2D array-like, shape=(n_frames, 3) order : int, default 2 tstep : float, default 1. tcorr : float, default 10000. norm : bool, default False |
---|
pytraj.
principal_axes
(traj=None, mask='*', dorotation=False, mass=True, top=None)¶compute principal axes
Parameters: | traj : Trajectory-like mask : str, default ‘*’ (all atoms) mass: bool, defaul True if `dorotation`, the system will be aligned along principal axes (apply for mutable system) top : Topology, optional |
---|---|
Returns: | out_0 : numpy array, shape=(n_frames, 3, 3) out_1: numpy array with shape=(n_frames, 3) |
pytraj.
randomize_ions
(traj, mask, around, by, overlap, seed=1, top=None, frame_indices=None)¶randomize_ions for given Frame with Topology
Parameters: | traj : Trajectory-like or a Frame
mask : str
frame_indices : {None, array-like}, optional top : Topology, optional (only needed if |
---|
Examples
>>> pt.randomize_ions(traj, mask='@Na+', around=':1-16', by=5.0, overlap=3.0, seed=113698)
pytraj.
crank
(data0, data1, mode='distance', dtype='ndarray')¶Parameters: | data0 : array-like data1 : array-like mode : str, {‘distance’, ‘angle’} dtype : str |
---|
Notes
Same as crank in cpptraj
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> distances = pt.distance(traj, [':3 :7', ':8 :12'])
>>> out = pt.crank(distances[0], distances[1])
pytraj.
autoimage
(traj, mask='', frame_indices=None, top=None)¶perform autoimage and return the updated-coordinate traj
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()[:]
>>> traj = pt.autoimage(traj)
pytraj.
closest
(traj=None, mask='*', solvent_mask=None, n_solvents=10, frame_indices=None, dtype='iterator', top=None)¶return either a new Trajectory or a frame iterator. Keep only n_solvents
closest to mask
Parameters: | traj : Trajectory-like | list of Trajectory-like/frames | frame iterator | chunk iterator mask: str, default ‘*’ (all solute atoms) top : Topology-like object, default=None, optional dtype : {‘iterator’, ‘trajectory’}, default ‘iterator’
|
---|---|
Returns: | out : (check above) |
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # obtain new traj, keeping only closest 100 waters
>>> # to residues 1 to 13 (index starts from 1) by distance to the first atom of water
>>> t = pt.closest(traj, mask='@CA', n_solvents=10)
pytraj.
search_neighbors
(traj=None, mask='', frame_indices=None, dtype='dataset', top=None)¶search neighbors
Returns: | pytraj.DatasetList, is a list of atom index arrays for each frame. Those arrays might not have the same lenghth |
---|
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> indices = pt.search_neighbors(traj, ':5<@5.0') # around residue 5 with 5.0 cutoff
pytraj.
watershell
(traj=None, solute_mask='', solvent_mask=':WAT', lower=3.4, upper=5.0, image=True, dtype='dataset', frame_indices=None, top=None)¶(adapted from cpptraj doc): Calculate numbers of waters in 1st and 2nd solvation shells (defined by <lower cut> (default 3.4 Ang.) and <upper cut> (default 5.0 Ang.)
Parameters: | traj : Trajectory-like solute_mask: solute mask solvent_mask: solvent mask lower : double, default 3.4
upper : double, default 5.0
image : bool, defaul True
dtype : return type, defaul ‘dataset’ top : Topology, optional |
---|
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.watershell(traj, solute_mask='!:WAT')
>>> data = pt.watershell(traj, solute_mask='!:WAT', lower=5.0, upper=10.)