pytraj.all_actions

note

  • just need to import pytraj and use the method
import pytraj as pt
pt.pca(traj)
  • Most of the methods’ names start with ‘calc_‘ but use can ignore ‘calc_‘. For example:

    pt.calc_pca(...) is the same as pt.pca(...)
    
pytraj.all_actions.translate(traj=None, command='', frame_indices=None, top=None)

translate coordinate

Examples

>>> import pytraj as pt
>>> # load to mutable trajectory by `load` method
>>> from pytraj.testing import get_fn
>>> fn, tn = get_fn('tz2')
>>> traj = pt.load(fn, tn)
>>> # do transform traj and return itself
>>> traj = pt.translate(traj, '@CA x 120.')
pytraj.all_actions.rotate(traj=None, command='', frame_indices=None, top=None)

Notes

rotate is an alias of do_rotation

Examples

>>> import pytraj as pt
>>> # load to mutable trajectory by `load` method
>>> from pytraj.testing import get_fn
>>> fn, tn = get_fn('tz2')
>>> traj = pt.load(fn, tn)
>>> # do transform traj and return itself
>>> traj = pt.rotate(traj, 'x 90')
pytraj.all_actions.autoimage(traj, mask='', frame_indices=None, top=None)

perform autoimage and return the coordinate-updated traj

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()[:]
>>> traj = pt.autoimage(traj)
pytraj.all_actions.scale(traj=None, command='', frame_indices=None, top=None)

Examples

>>> import pytraj as pt
>>> # load to mutable trajectory by `load` method
>>> from pytraj.testing import get_fn
>>> fn, tn = get_fn('tz2')
>>> traj = pt.load(fn, tn)
>>> # do transform traj and return itself
>>> traj = pt.scale(traj, '@CA x 1.2')
pytraj.all_actions.calc_distance(traj=None, mask='', frame_indices=None, dtype='ndarray', top=None, image=False, n_frames=None)

calculate 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

only need to provide n_frames if traj does not have this info

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.all_actions.calc_dihedral(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None, *args, **kwargs)

calculate dihedral 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 dihedral angle for four atoms, using amber mask
>>> pt.dihedral(traj, '@1 @3 @10 @20')
array([ 23.32244362,  23.5386922 ,  14.26831569,  14.58865946,
        23.98675475,  26.18419185,   6.06982926,  13.57158505,
        16.59013076,  29.99131573])
>>> # calculate dihedral angle for four groups of atoms, using amber mask
>>> data = pt.dihedral(traj, '@1,37,8 @2,4,6 @5,20 @21,22')
>>> # calculate dihedral angle for four residues, using amber mask
>>> data = pt.dihedral(traj, ':1 :10 :20 :22')
>>> # calculate multiple dihedral angles for four residues, using amber mask
>>> # angle for residue 1, 10, 20, 30; angle between residue 3, 20, 30, 40
>>> # (when using atom string mask, index starts from 1)
>>> pt.dihedral(traj, [':1 :10 :20 :30', ':3 :20 :30 :40'])
array([[-166.81829116, -160.19029669, -158.56560062, ..., -169.10064772,
        -158.81655586, -165.28873555],
       [  -0.60994639,    0.78261235,    1.86394369, ...,    1.21170489,
          -1.16762607,   -3.08412049]])
>>> # calculate dihedral angle for a series of atoms, using array for atom mask
>>> # dihedral angle for atom 1, 5, 8, 10, dihedral for atom 4, 10, 20, 30 (index starts from 0)
>>> data = pt.dihedral(traj, [[1, 5, 8, 10], [4, 10, 20, 30]])
pytraj.all_actions.calc_radgyr(traj=None, mask='', top=None, nomax=True, frame_indices=None, dtype='ndarray')

calc radgyr

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.all_actions.calc_angle(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None, *args, **kwargs)

calculate 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.all_actions.calc_surf(traj=None, mask='', dtype='ndarray', frame_indices=None, top=None)

calc surf (LCPO method)

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.surf(traj, '@CA')
pytraj.all_actions.calc_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.all_actions.calc_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.calc_volume(traj, '@CA')
pytraj.all_actions.calc_matrix(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None)

compute different type of matrices

Parameters:

traj : Trajectory-like

mask : str, type of matrix and atom mask

top : Topology, optional

dtype: return data type

frame_indices : {None, array-like}

if not None, perform calculation for given frame indices

Notes

If user wants to use specify matrix’s method name, see also pytraj.matrix

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_trpcage()
>>> mat = pt.calc_matrix(traj, 'covar @CA')
>>> # this is equal to
>>> mat2 = pt.matrix.covar(traj, '@CA')
>>> import numpy as np
>>> np.testing.assert_equal(mat, mat2)
pytraj.all_actions.calc_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 “”

cpptraj’s command/mask

kfile : str, default None, optional

Dir for Karplus file. If “None”, use $AMBERHOME dir

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.all_actions.calc_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

lower cut distance

upper : double, default 5.0

upper cut distance

image : bool, defaul True

do autoimage if 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.)
pytraj.all_actions.calc_vector(traj=None, command='', frame_indices=None, dtype='ndarray', top=None)

perform vector calculation. See example below. Same as ‘vector’ command in cpptraj.

Parameters:

traj : Trajectory-like or iterable that produces pytraj.Frame

command : str or a list of strings, cpptraj command

frame_indices : array-like, optional, default None

only perform calculation for given frame indices

dtype : output’s dtype, default ‘ndarray’

top : Topology, optional, default None

Returns:

out : numpy ndarray, shape (n_frames, 3) if command is a string

numpy ndarray, shape (n_vectors, n_frames, 3) if command is a list of strings

Notes

It’s faster to calculate with a list of commands. For example, if you need to perform 3 calculations for ‘ucellx’, ‘boxcenter’, ‘box’ like below:

>>> data = pt.calc_vector(traj, "ucellx")
>>> data = pt.calc_vector(traj, "boxcenter")
>>> data = pt.calc_vector(traj, "box")

You should use a list of commands for faster calculation.

>>> comlist = ['ucellx', 'boxcenter', 'box']
>>> data = pt.calc_vector(traj, comlist)

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.calc_vector(traj, "@CA @CB")
>>> data = pt.calc_vector(traj, [("@CA @CB"),])
>>> data = pt.calc_vector(traj, "principal z")
>>> data = pt.calc_vector(traj, "principal x")
>>> data = pt.calc_vector(traj, "ucellx")
>>> data = pt.calc_vector(traj, "boxcenter")
>>> data = pt.calc_vector(traj, "box")
pytraj.all_actions.calc_multivector(traj, resrange, names, top=None, dtype='dataset', frame_indices=None)
Parameters:

traj : Trajectory-like

resrange : str, residue range

names : {str, tuple of str}

top : Topology, optional

dtype : str, default ‘dataset’

frame_indices : {None, 1D array-like}, optional, default None

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> vecs = pt.calc_multivector(traj, resrange='1-5', names=('C', 'N'))
>>> vecs = pt.calc_multivector(traj, resrange='1-5', names='C N')
pytraj.all_actions.calc_volmap(traj, mask, grid_spacing, size=None, center=None, buffer=3.0, centermask='*', radscale=1.36, peakcut=0.05, top=None, dtype='ndarray', frame_indices=None)

(combined with cpptraj doc) Grid data as a volumetric map, similar to the volmap command in VMD. The density is calculated by treating each atom as a 3-dimensional Gaussian function whose standard deviation is equal to the van der Waals radius

Parameters:

mask : {str, array-like}, default all atoms

the atom selection from which to calculate the number density

grid_spacing : tuple, grid spacing in X-, Y-, Z-dimensions, require

size : {None, tuple}, default None

if tuple, size must have length of 3

center : {None, tuple}, default None

if not None, center is tuple of (x, y, z) of center point

buffer : float, default 3.0 Angstrom

buffer distance (Angstrom), by which the edges of the grid should clear every atom of the centermask (or default mask if centermask is omitted) in every direction. The buffer is ignored if the center and size are specified.

centermask : str

radscale : float, default 1.36 (to match to VMD calculation)

factor by which to scale radii (by devision)

peakcut : float

Examples

>>> import pytraj as pt
>>> # load all frames to memory
>>> traj = pt.datafiles.load_tz2_ortho()[:]
>>> # do fitting and centering before perform volmap
>>> traj = traj.superpose(mask=':1-13').center(':1-13 mass origin')
>>> data = pt.volmap(traj, mask=':WAT@O', grid_spacing=(0.5, 0.5, 0.5), buffer=2.0, centermask='!:1-13', radscale=1.36)
pytraj.all_actions.calc_rdf(traj=None, solvent_mask=':WAT@O', solute_mask='', maximum=10.0, bin_spacing=0.5, image=True, density=0.033456, volume=False, center_solvent=False, center_solute=False, intramol=True, frame_indices=None, top=None)

compute radial distribtion function. Doc was adapted lightly from cpptraj doc

Parameters:

traj : Trajectory-like, require

solvent_mask : solvent mask, default None, required

bin_spacing : float, default 0.5, optional

bin spacing

maximum : float, default 10., optional

max bin value

solute_mask: str, default None, optional

if specified, calculate RDF of all atoms in solvent_mask to each atom in solute_mask

image : bool, default True, optional

if False, do not image distance

density : float, default 0.033456 molecules / A^3, optional

volume : determine density for normalization from average volume of input frames

center_solvent : bool, default False, optional

if True, calculate RDF from geometric center of atoms in solvent_mask to all atoms in solute_mask

center_solute : bool, default False, optional

if True, calculate RDF from geometric center of atoms in solute_mask to all atoms in solvent_mask

intramol : bool, default True, optional

if False, ignore intra-molecular distances

frame_indices : array-like, default None, optional

Returns:

a tuple of bin_centers, rdf values

Notes

  • install pytraj and libcpptraj with openmp to speed up calculation
  • do not use this method with pytraj.pmap

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.rdf(traj, solvent_mask=':WAT@O', bin_spacing=0.5, maximum=10.0, solute_mask=':WAT@O')
>>> data[0]
array([ 0.25,  0.75,  1.25, ...,  8.75,  9.25,  9.75])
>>> data[1]
array([ 0.        ,  0.        ,  0.        , ...,  0.95620052,
        0.95267934,  0.95135242])
>>> # use array-like mask
>>> atom_indices = pt.select(':WAT@O', traj.top)
>>> data = pt.rdf(traj, solvent_mask=':WAT@O', bin_spacing=0.5, maximum=10.0, solute_mask=atom_indices)
pytraj.all_actions.calc_pairdist(traj, mask='*', mask2='', delta=0.1, dtype='ndarray', top=None, frame_indices=None)

compute pair distribution function

Parameters:

traj : Trajectory-like

mask : str, default all atoms

delta : float, default 0.1

dtype : str, default ‘ndarray’

dtype of return data

top : Topology, optional

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.calc_pairdist(traj)
pytraj.all_actions.calc_multidihedral(traj=None, dhtypes=None, resrange=None, define_new_type=None, range360=False, dtype='dataset', top=None, frame_indices=None)

perform dihedral search

Parameters:

traj : Trajectory-like object

dhtypes : dihedral type, default None

if None, calculate all supported dihedrals

resrange : str | array-like

residue range for searching. If resrange is string, use index starting with 1 (cpptraj convertion) if resrange is array-like, use index starting from 0 (python convention)

define_new_type : str

define new type for searching

range360 : bool, default False

if True: use 0-360

top : Topology | str, optional

only need to have ‘top’ if can not find it in traj

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, dhtypes='phi psi')
>>> data = pt.multidihedral(traj, dhtypes='phi psi', resrange='3-7')
>>> data = pt.multidihedral(traj, dhtypes='phi psi', resrange=[3, 4, 8])
pytraj.all_actions.calc_atomicfluct(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None)

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.atomicfluct(traj, '@CA')
>>> data[:3]
array([[  5.        ,   0.61822273],
       [ 16.        ,   0.5627449 ],
       [ 40.        ,   0.53717119]])
pytraj.all_actions.calc_COM(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.
>>> pt.calc_center_of_mass(traj(stop=2), ':3')
array([[-0.661702  ,  6.69124347,  3.35159413],
       [ 0.5620708 ,  7.82263042, -0.72707798]])
pytraj.all_actions.calc_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.
>>> pt.calc_center_of_mass(traj(stop=2), ':3')
array([[-0.661702  ,  6.69124347,  3.35159413],
       [ 0.5620708 ,  7.82263042, -0.72707798]])
pytraj.all_actions.calc_center_of_geometry(traj=None, mask='', top=None, dtype='ndarray', frame_indices=None)
pytraj.all_actions.calc_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

if mask is “”, use all atoms

metric : {‘rms’, ‘dme’, ‘srmsd’, ‘nofit’}

if ‘rms’, perform rms fit if ‘dme’, use distance RMSD if ‘srmsd’, use symmetry-corrected RMSD if ‘nofit’, perform rmsd without fitting

top : Topology, optional, default=None

dtype: ndarray

return type

mat_type : str, {‘full’, ‘half’}

if ‘full’: return 2D array, shape=(n_frames, n_frames) if ‘half’: return 1D array, shape=(n_frames*(n_frames-1)/2, )

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.all_actions.calc_grid(traj=None, command='', top=None, dtype='dataset')
pytraj.all_actions.calc_atomiccorr(traj, mask='', cut=None, min_spacing=None, byres=True, frame_indices=None, dtype='ndarray', top=None)

compute average correlations between the motion of atoms in mask.

Parameters:

traj : Trajectory-like

mask : atom mask

cut : {None, float}, default None

if not None, only print correlations with absolute value greater than cut

min_spacing : {None, float}, default None

if not None, only calculate correlations for motion vectors spaced min_spacing apart

byres : bool, default False

if False, compute atomic motion vetor if True, Calculate motion vectors for entire residues (selected atoms in residues only).

pytraj.all_actions.calc_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])

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.calc_bfactors(traj, byres=True)
pytraj.all_actions.calc_diffusion(traj, mask='', tstep=1.0, individual=False, top=None, dtype='dataset', frame_indices=None)

compute diffusion

Parameters:

traj : Trajectory-like or iterator

mask : str, default ‘’ (all atoms)

tstep : float, time step between frames, default 1.0 ps

individual : bool, default False

top : Topology, optional

dtype : str, default ‘dataset’

frame_indices : array or None

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.diffusion(traj, dtype='dict')
>>> data['X']
array([ 0.        ,  0.87027302,  1.64626022,  2.26262651,  2.98068114,
        3.57075535,  4.07030655,  4.71894512,  5.42302306,  6.01310377])
pytraj.all_actions.calc_distance_rmsd(traj=None, mask='', ref=0, top=None, dtype='ndarray', frame_indices=None)

compute distance rmsd between traj and reference

Parameters:

traj : Trajectory-like or iterator that produces Frame

ref : {int, Frame}, default 0 (1st Frame)

mask : str

top : Topology or str, optional, default None

dtype : return dtype, default ‘ndarray’

Returns:

1D ndarray if dtype is ‘ndarray’ (default)

Examples

>>> import pytraj as pt
>>> # compute distance_rmsd to last frame
>>> traj = pt.datafiles.load_tz2_ortho()
>>> data = pt.distance_rmsd(traj, ref=-1)
>>> # compute distance_rmsd to first frame with mask = '@CA'
>>> data = pt.distance_rmsd(traj, ref=0, mask='@CA')
pytraj.all_actions.calc_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.all_actions.calc_pairwise_distance(traj=None, mask_1='', mask_2='', top=None, dtype='ndarray', frame_indices=None)

compute pairwise distance between atoms in mask_1 and atoms in mask_2

Parameters:

traj : Trajectory-like or iterable that produces Frame

mask_1: string or 1D array-like

mask_2: string or 1D array-like

...

Returns:

out_1 : numpy array, shape (n_frames, n_atom_1, n_atom_2)

out_2 : atom pairs, shape=(n_atom_1, n_atom_2, 2)

Notes

This method is only fast for small number of atoms.

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> mat = pt.calc_pairwise_distance(traj, '@CA', '@CB')
pytraj.all_actions.calc_rmsd_nofit(traj=None, mask='', ref=0, mass=False, frame_indices=None, top=None, dtype='ndarray', **kwd)

See also

calc_rmsd

pytraj.all_actions.calc_rotation_matrix(traj=None, mask='', ref=0, mass=False, frame_indices=None, top=None, with_rmsd=False)

compute rotation matrix with/without rmsd

Parameters:

traj : Trajectory-like

ref : {int, Frame}, default 0 (first Frame)

reference

mask : str, default all atoms

mass : bool, default False

if True, rmsfit with mass

frame_indices : {None, array-like}

if not None, compute for given indices

top : Topology, optional

with_rmsd : bool, default False

  • if False, return only rotation matrix.
  • if True, return rotation matrix and rmsd values
Returns:

out : if with_rmsd=False, return numpy array, shape (n_frames, 3, 3)

if with_rmsd=True, return a tuple (mat, rmsd)

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> mat = pt.calc_rotation_matrix(traj, mask='@CA')
>>> mat.shape
(101, 3, 3)
pytraj.all_actions.calc_pca(traj, mask, n_vecs=2, fit=True, ref=None, ref_mask=None, dtype='ndarray', top=None)

perform PCA analysis by following below steps:

  • (optional) perform rmsfit to reference if needed
  • compute covariance matrix
  • diagonalize the matrix to get eigenvectors and eigenvalues
  • perform projection of each frame with mask to each eigenvector
Parameters:

traj : Trajectory

traj must be pytraj.Trajectory, which can be created by pytraj.load method.

mask : str

atom mask for covariance matrix and projection

n_vecs : int, default 2

number of eigenvectors. If user want to compute projection for all eigenvectors, should specify n_vecs=-1 (or a negative number)

fit : bool, default True

if True, perform fitting before compute covariance matrix if False, no fitting (keep rotation and translation). In this case, pytraj will ignore ref argument.

ref : {None, Frame, int}, default None

if None, trajectory will be superposed to average structure if is Frame or integer value, trajectory will be superposed to given reference

ref_mask : {None, str}, default None (use mask)

if None, use mask for fitting if str, use this given mask for fitting

dtype : return datatype

top : Topology, optional

Returns:

out1: projection_data, ndarray with shape=(n_vecs, n_frames)

out2: tuple of (eigenvalues, eigenvectors)

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()[:]
>>> # compute pca for first and second modes
>>> pca_data = pt.pca(traj, '!@H=', n_vecs=2)
>>> # get projection values
>>> pca_data[0] 
array([[  4.93425131,  13.80002308,  20.61605835, ..., -57.92280579,
        -61.25728607, -52.85142136],
       [  4.03333616,  -6.9132452 , -14.53991318, ...,  -6.757936  ,
          2.1086719 ,  -3.60922861]], dtype=float32)
>>> # get eigenvalues for first 2 modes
>>> pca_data[1][0] 
array([ 1399.36472919,   240.42342439])
>>> # compute pca for all modes
>>> pca_data = pt.pca(traj, '!@H=', n_vecs=-1)
>>> # does not perform fitting
>>> data = pt.pca(traj, mask='!@H=', fit=False)
>>> # provide different mask for fitting
>>> data = pt.pca(traj, mask='!@H=', fit=True, ref=0, ref_mask='@CA')
pytraj.all_actions.get_average_frame(traj, mask='', frame_indices=None, dtype='frame', autoimage=False, rmsfit=None, top=None)

get mean structure for a given mask and given frame_indices

Parameters:

traj : Trajectory-like or iterable that produces Frame

mask : None or str, default None (all atoms)

frame_indices : iterable that produces integer, default None, optional

frame indices

dtype: str, {‘frame’, ‘traj’}, default ‘frame’

return type, either Frame (does not have Topology information) or ‘traj’

autoimage : bool, default False

if True, performa autoimage

rmsfit : object, {Frame, int, tuple, None}, default None

if rmsfit is not None, perform rms fit to reference.

Notes

if autoimage=True and having rmsfit, perform autoimage first and do rmsfit

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # get average structure from whole traj, all atoms
>>> frame = pt.mean_structure(traj)
>>> # get average structure from several frames, '@CA' atoms"
>>> frame = pt.mean_structure(traj, '@CA', frame_indices=range(2, 8, 2))
>>> # get average structure but do autoimage and rmsfit to 1st Frame
>>> frame = pt.mean_structure(traj(autoimage=True, rmsfit=0))
>>> # get average structure but do autoimage and rmsfit to 1st Frame.
>>> frame = pt.mean_structure(traj(autoimage=True, rmsfit=0, frame_indices=[0, 5, 6]))
pytraj.all_actions.get_velocity(traj, mask=None, frame_indices=None)

get velocity for specify frames with given mask

Parameters:

traj : Trajectory-like or iterable that produces Frame

mask : str, default None (use all atoms), optional

atom mask

frame_indices : iterable that produces integer, default None, optional

if not None, only get velocity for given frame indices

Notes

Since pytraj has limited support for force and velocity info, if user wants to load both from disk, should iterate the TrajectoryIterator (got by pytraj.iterload method)

import pytraj as pt
forces = []
velocities = []

traj = pt.iterload(filename, topology_filename)

for frame in traj:
    forces.append(frame.force)
    velocities.append(frame.velocity)

# Note: pytraj can efficiently load arbitary frame indices to memory
for frame in pt.iterframe(traj, frame_indices=[0, 8, 8, 100, 1000]): pass

Examples

>>> vels = pt.get_velocity(traj, frame_indices=[0, 3]) 
pytraj.all_actions.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

Atom mask/indices

ref : {Frame, int}, default=0 (first frame)

Reference frame or index.

ref_mask: str, optional

if given, use it instead of mask

nofit : bool, default False

if False, perform fitting (rotation and translation). if traj is mutable, its coordinates will be updated if True, not fitting.

mass : bool, default False

if True, include mass

update_coordinate : bool, default True

if True, coordinates will be updated. But this only apply to mutable Trajectory if False (same as nomod in cpptraj), no modification

frame_indices : int 1D array-like, default None

if not None, only compute rmsd for given frame indices

top : {Topology, str}, default None, optional

dtype : return data type, default=’ndarray’

Notes

if traj is mutable, its coordinates will be updated.

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])
pytraj.all_actions.align_principal_axis(traj=None, mask='*', top=None, frame_indices=None, mass=False)

Notes

apply for mutatble traj (Trajectory, Frame)

pytraj.all_actions.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.all_actions.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’

if ‘iterator’, return a tuple of Frame iterator and new Toplogy. ‘iterator’ is good for streaming large trajectory. if ‘trajectory’, return a new Trajectory.

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.all_actions.transform(traj, by, frame_indices=None)

transform pytraj.Trajectory by a series of cpptraj’s commands

Parameters:

traj : Mutable Trajectory

by : list of cpptraj commands

frame_indices : {None, array-like}, default None

if not None, perform tranformation for specific frames.

Returns:

transformed Trajectory. Trajectory’s coordinates will be inplace-updated

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # perform 'autoimage', then 'rms', then 'center'
>>> traj = pt.transform(traj[:], by=['autoimage', 'rms', 'center :1-5'])
pytraj.all_actions.native_contacts(traj=None, mask='', mask2='', ref=0, dtype='dataset', distance=7.0, image=True, include_solvent=False, byres=False, frame_indices=None, top=None)

compute native contacts

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # use 1st frame as reference, don't need specify ref Frame
>>> data = pt.native_contacts(traj)
>>> # explicitly specify reference, specify distance cutoff
>>> ref = traj[3]
>>> data = pt.native_contacts(traj, ref=ref, distance=8.0)
>>> # use integer array for mask
>>> data = pt.native_contacts(traj, mask=range(100), mask2=[200, 201], ref=ref, distance=8.0)
pytraj.all_actions.set_dihedral(traj, resid='1', dihedral_type=None, deg=0, top=None)
Returns:updated traj

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> # make mutable traj by loading all frames to disk
>>> traj = traj[:]
>>> traj = pt.set_dihedral(traj, resid='3', dihedral_type='phi', deg=60)
>>> traj = pt.set_dihedral(traj, resid=2, dihedral_type='phi', deg=60)
pytraj.all_actions.check_structure(traj, mask='', options='', frame_indices=None, top=None, dtype='ndarray')

check if the structure is ok or not

Parameters:

traj : Trajectory-like

mask: str, default all atoms

options : str, default ‘’

extra cpptraj options

dtype : str, default ‘ndarray’

Returns:

out : 1D-array

number of problems for each frame

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_rna()
>>> failures = pt.check_structure(traj[0], top=traj.top)
pytraj.all_actions.mean_structure(traj, mask='', frame_indices=None, dtype='frame', autoimage=False, rmsfit=None, top=None)

get mean structure for a given mask and given frame_indices

Parameters:

traj : Trajectory-like or iterable that produces Frame

mask : None or str, default None (all atoms)

frame_indices : iterable that produces integer, default None, optional

frame indices

dtype: str, {‘frame’, ‘traj’}, default ‘frame’

return type, either Frame (does not have Topology information) or ‘traj’

autoimage : bool, default False

if True, performa autoimage

rmsfit : object, {Frame, int, tuple, None}, default None

if rmsfit is not None, perform rms fit to reference.

Notes

if autoimage=True and having rmsfit, perform autoimage first and do rmsfit

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # get average structure from whole traj, all atoms
>>> frame = pt.mean_structure(traj)
>>> # get average structure from several frames, '@CA' atoms"
>>> frame = pt.mean_structure(traj, '@CA', frame_indices=range(2, 8, 2))
>>> # get average structure but do autoimage and rmsfit to 1st Frame
>>> frame = pt.mean_structure(traj(autoimage=True, rmsfit=0))
>>> # get average structure but do autoimage and rmsfit to 1st Frame.
>>> frame = pt.mean_structure(traj(autoimage=True, rmsfit=0, frame_indices=[0, 5, 6]))
pytraj.all_actions.lowestcurve(data, points=10, step=0.2)

compute lowest curve for data

pytraj.all_actions.make_structure(traj=None, mask='', top=None)
pytraj.all_actions.replicate_cell(traj=None, mask='', direction='all', top=None)

create a trajectory where the unit cell is replicated in 1 or more direction (up to 27)

Parameters:

traj : Trajectory-like or Frame iterator

mask : str, default: “”

if default, using all atoms else: given mask

direction: {‘all’, ‘dir’} or list/tuple of <XYZ> (below)

if ‘all’, replicate cell once in all possible directions if ‘dir’, need to specify the direction with format ‘dir <XYZ>’, where each X (Y, Z) is either 0, 1 or -1 (see example below)

top : Topology, optional, default: None

Returns:

traj : pytraj.Trajectory

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> new_traj = pt.replicate_cell(traj, direction='all')
>>> new_traj = pt.replicate_cell(traj, direction='dir 001 dir 111')
>>> new_traj = pt.replicate_cell(traj, direction='dir 001 dir 1-10')
>>> new_traj = pt.replicate_cell(traj, direction='dir 001 dir 1-10')
>>> # similiar usage
>>> new_traj = pt.replicate_cell(traj, direction=('001', '0-10'))
>>> new_traj = pt.replicate_cell(traj, direction=['001', '0-10'])
pytraj.all_actions.pucker(traj=None, pucker_mask=("C1'", "C2'", "C3'", "C4'", "O4'"), resrange=None, top=None, dtype='dataset', range360=False, method='altona', use_com=True, amplitude=False, offset=None)

compute pucker

Parameters:

traj : Trajectory-like

pucker_mask : str

resrange : None or array of int

top : Topology, optional

dtype : str, return type

range360: bool, use 360 or 180 scale

method : {‘altona’, ‘cremer’}, default ‘altona’

use_com : bool

amplitude : bool, default False

offset : None or float

Returns:

Dataset

pytraj.all_actions.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)

out[0]: regular rmsd out[1:]: perres rmsd for all given residues out.values will return corresponding numpy array

pytraj.all_actions.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

traj must be mutable

mask : str

cpptraj command

frame_indices : {None, array-like}, optional

top : Topology, optional (only needed if traj does not have it)

Examples

>>> pt.randomize_ions(traj, mask='@Na+', around=':1-16', by=5.0, overlap=3.0, seed=113698) 
pytraj.all_actions.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.all_actions.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.all_actions.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.all_actions.acorr(data, dtype='ndarray', option='')

compute autocorrelation

Parameters:

data : 1d array-like

dtype: return type, default ‘ndarray’

covar : bool, default True

option : str

more cpptraj options

Notes

Same as autocorr in cpptraj

pytraj.all_actions.projection(traj, mask='', eigenvalues=None, eigenvectors=None, scalar_type='covar', average_coords=None, frame_indices=None, dtype='ndarray', top=None)

compute projection along given eigenvectors

Parameters:

traj : Trajectory-like

mask : atom mask, either string or array-like

eigenvalues : 1D array-like

eigenvectors : 2D array-like

scalar_type : str, {‘covar’, ‘mwcovar’, }, default ‘covar’

make sure to provide correct scalar_type. Note: not yet support ‘dihcovar’ and ‘idea’

average_coords : 3D array-like, optional, default None

average coordinates. If ‘None’, pytraj will compute mean_structure with given mask

frame_indices : array-like

If not None, compute projection for given frames.

dtype : str, return data type, default ‘ndarray’

top : Topology, optional, default None

Returns:

projection_data : ndarray, shape=(n_vecs, n_frames)

Examples

>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()
>>> mat = pt.matrix.covar(traj, '@CA')
>>> eigenvalues, eigenvectors = pt.matrix.diagonalize(mat, 2)
>>> # since we compute covariance matrix, we need to specify
>>> # scalar_type = 'covar'
>>> scalar_type = 'covar'
>>> data = pt.projection(traj, '@CA', eigenvalues=eigenvalues, eigenvectors=eigenvectors, scalar_type=scalar_type)
>>> data.shape
(2, 101)