Note
Just need to import pytraj and use the method directly
import pytraj as pt
pt.pca(traj)
pt.surf(traj, mask='@CA')
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.
analyze_modes
(mode_type, eigenvectors, eigenvalues, scalar_type='mwcovar', options='', dtype='dict')¶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 updated-coordinate traj
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()[:]
>>> traj = pt.autoimage(traj)
pytraj.all_actions.
image
(traj, mask='', frame_indices=None, top=None)¶perform imaging and return the updated-coordinate traj
Notes
User should always try to use autoimage first
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()[:]
>>> traj = pt.image(traj, 'origin center :WAT')
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.
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
dtype: str, {‘frame’, ‘traj’}, default ‘frame’
autoimage : bool, default False
rmsfit : object, {Frame, int, tuple, None}, default None
|
---|
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
frame_indices : iterable that produces integer, default None, optional
|
---|---|
Returns: | out : 3D numpy array, shape (n_frames, n_atoms, 3) |
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.
atom_map
(traj, ref, rmsfit=False)¶Limited support for cpptraj atommap
Parameters: | traj : Trajectory-like ref : Trajectory-like with one frame rmsfit : bool, default False
|
---|---|
Returns: | out : Tuple[str, np.ndarray]
|
Notes
This method in pytraj is not mature yet.
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’
|
---|---|
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
|
---|---|
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=0, 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=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 ‘’
dtype : str, default ‘ndarray’ |
---|---|
Returns: | out : Tuple[np.ndarray, str]
|
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_rna()
>>> failures = pt.check_structure(traj[:1])
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
dtype: str, {‘frame’, ‘traj’}, default ‘frame’
autoimage : bool, default False
rmsfit : object, {Frame, int, tuple, None}, default None
|
---|
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
Parameters: | data : 2D array-like points : number of lowest points in each bin, default 10 step : step size, default 0.2 |
---|---|
Returns: | 2d array |
pytraj.all_actions.
make_structure
(traj, command='', ref=None)¶limited support for make_structure
Parameters: | traj : Trajectory command : cpptraj command or a list of cpptraj commands ref : None or a Frame or a Trajectory with one Frame
|
---|---|
Returns: | traj : itself |
Notes
Apply dihedrals to specified residues using arguments found in <List of Args>, where an argument is 1 or more of the following arg types:
‘<sstype>:<res range>’
Apply SS type (phi/psi) to residue range.
<sstype> standard = alpha, left, pp2, hairpin, extended
- <sstype> turn = typeI, typeII, typeVIII, typeI’, typeII, typeVIa1, typeVIa2, typeVIb
- Turns are applied to 2 residues at a time, so resrange must be divisible by 4.
‘<custom ss>:<res range>:<phi>:<psi>’
Apply custom <phi>/<psi> to residue range.
‘<custom turn>:<res range>:<phi1>:<psi1>:<phi2>:<psi2>’
Apply custom turn <phi>/<psi> pair to residue range.
‘<custom dih>:<res range>:<dih type>:<angle>’
Apply <angle> to dihedrals in range.
<dih type> = alpha beta gamma delta epsilon zeta nu1 nu2 h1p c2p chin phi psi chip omega
‘<custom dih>:<res range>:<at0>:<at1>:<at2>:<at3>:<angle>[:<offset>]’
Apply <angle> to dihedral defined by atoms <at1>, <at2>, <at3>, and <at4>.
- Offset -2=<a0><a1> in previous res, -1=<a0> in previous res,
- 0=All <aX> in single res, 1=<a3> in next res, 2=<a2><a3> in next res.
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2()[:]
>>> # Make alpha helix for residue 1 to 12)
>>> traj = pt.make_structure(traj, "alpha:1-12")
>>> # Make hairpin for residue 1 to 5 and make alpha helix for residue 6 to 12
>>> traj = pt.make_structure(traj, ["hairpin:1-5", "alpha:6-12"])
>>> # From cpptraj example:
>>> # Make residues 1 and 12 'extended', residues 6 and 7 a type I' turn, and two
>>> # custom assignments, one (custom1) for residues 2-5, the other (custom2) for residues 8-11:
>>> traj = pt.make_structure(traj, ["extended:1,12",
... "custom1:2-5:-80.:130.:-130.:140.",
... "typeI':6-7",
... "custom2:8-11:-140.:170.:-100.:140."])
>>>
>>> # Make new structure from reference
>>> def make_new_by_using_ref():
... tz2_parm7 = 'tz2.parm7'
... ref_rst7 = 'tz2.rst7'
... trajin_rst7 = pp2.rst7.save'
... traj = pt.load(trajin_rst7, top=tz2_parm7)
... ref = pt.load(ref_rst7, top=tz2_parm7)
... # make dihedral angle for targeted residues 1 to 13 by using
... # reference residues 1 to 13
... pt.make_structure(traj, "ref:1-13:1-13", ref=ref)
pytraj.all_actions.
replicate_cell
(traj=None, mask='', direction='all', frame_indices=None, 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: “”
direction: {‘all’, ‘dir’} or list/tuple of <XYZ> (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)
|
---|
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
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.all_actions.
velocityautocorr
(traj, mask='', maxlag=-1, tstep=1.0, direct=True, norm=False, usevelocity=False, dtype='ndarray', top=None, velocity_arr=None)¶Parameters: | traj : Trajectory-like mask : str
maxlag : int, default -1
tstep : float, default 1.0 direct : bool, default True
norm : bool, default False
usevelocity : bool, default False
dtype : str, default ‘ndarray’
top : None or Topology, default None, optional velocity_arr : None or 3D like array, default None
|
---|
Notes
If you create Trajectory by pytraj.load method, there is no velocity information. So if you want to use usevelocity=True, you need to provide 3D-array velocity_arr
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
|
---|
Notes
Same as autocorr in cpptraj
pytraj.all_actions.
projection
(traj, mask='', eigenvectors=None, eigenvalues=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’
average_coords : 3D array-like, optional, default None
frame_indices : array-like
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')
>>> eigenvectors, eigenvalues = 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)
pytraj.all_actions.
superpose
(traj, *args, **kwd)¶pytraj.all_actions.
strip
(obj, mask)¶return a new Trajectory or FrameIterator or Topology with given mask.
Notes
This method is trying to be smart. If you give it an in-memory Trajectory, it will return a corresponding in-memory Trajectory. If you give it an out-of-memory TrajectoryIterator, it will give you a corresponding FrameIterator object (out-of-memory). If you give it a Topology, it will return a new stripped Topology.
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> traj.n_atoms
5293
>>> pt.strip(traj, '!@CA').n_atoms
12
>>> pt.strip(traj.top, '!@CA').n_atoms
12
>>> pt.strip('!@CA', traj.top).n_atoms
12
>>> t0 = traj[:3]
>>> pt.strip(t0, '!@CA').n_atoms
12
>>> fi = pt.iterframe(traj, stop=3)
>>> fi = pt.strip(fi, '!@CA')
>>> for frame in fi:
... print(frame)
<Frame with 12 atoms>
<Frame with 12 atoms>
<Frame with 12 atoms>
>>> # raise ValueError
>>> pt.strip(0, '@CA')
Traceback (most recent call last):
...
ValueError: object must be either Trajectory or Topology
pytraj.all_actions.
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.all_actions.
gist
(traj, grid_center=[0.0, 0.0, 0], grid_dim=[40, 40, 40], grid_spacing=0.5, do_order=False, do_eij=False, reference_density=0.0334, temperature=300.0, options='', dtype='dict')¶minimal support for gist command in cpptraj
Parameters: | traj : Trajectory-like grid_center : 1-D array-like or str, default [0., 0., 0.] (origin)
grid_dim: 1-D array-like or str, default [40, 40, 40]
grid_spacing: float, default 0.5 do_order : bool, default False do_eij : bool, default False reference_density : float, default 0.0334
options : str
temperature : float, default 300. dtype : str, default ‘dict’
|
---|---|
Returns: | out : dict (or another data type based on dtype)
|
Notes
Syntax might be changed. There is a bug in pytraj that causes segmentation fault sometimes.
pytraj.all_actions.
center
(traj=None, mask='', center='box', mass=False, top=None, frame_indices=None)¶Center coordinates in mask to specified point.
Parameters: | traj : Trajectory-like or Frame iterator mask : str, mask center : str, {‘box’, ‘origin’, array-like}, default ‘box’
mass : bool, default: False
top : Topology, optional, default: None |
---|---|
Returns: | updated traj |
See also
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_tz2_ortho()
>>> # load all frames to memory so we can 'mutate' them
>>> traj = traj[:]
>>> # all atoms, center to box center (x/2, y/2, z/2)
>>> traj = pt.center(traj)
>>> # center at origin, use @CA
>>> traj = pt.center(traj, '@CA', center='origin')
>>> # center to box center, use mass weighted
>>> traj = pt.center(traj, mass=True)
>>> traj = pt.center(traj, ':1', mass=True)
pytraj.all_actions.
wavelet
(traj, command)¶wavelet analysis
Parameters: | traj : Trajectory-like command : str, cpptraj command |
---|---|
Returns: | out : dict |
Notes
you need to type cpptraj command. Please check cpptraj manual for further info if you really want to use it.
so this method is only good for small trajectory that fit to your RAM.
version added: 1.0.6
Examples
>>> import pytraj as pt
>>> traj = pt.datafiles.load_dpdp()
>>> c0 = 'nb 10 s0 2 ds 0.25 type morlet correction 1.01 chival 0.25 :1-22'
>>> c1 = 'cluster minpoints 66 epsilon 10.0'
>>> command = ' '.join((c0, c1))
>>> wavelet_dict = pt.wavelet(traj, command)
pytraj.all_actions.
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)
mask : str, default all atoms mass : bool, default False
frame_indices : {None, array-like}
top : Topology, optional with_rmsd : bool, default False
|
---|---|
Returns: | out : if with_rmsd=False, return numpy array, shape (n_frames, 3, 3)
|
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.
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.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)
|
---|
pytraj.all_actions.
rmsd_nofit
(traj=None, mask='', ref=0, mass=False, frame_indices=None, top=None, dtype='ndarray', **kwd)¶compute rmsd without fitting (translating and rotating)
Parameters: | traj : Trajectory-like mask : str ref : Frame or int mass : bool, default False
frame_indices : 1D array-like, default None
|
---|
Notes
This method is equal to pytraj.rmsd(traj, mask, ref, nofit=True, ...)
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
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.all_actions.
symmrmsd
(traj, mask='', ref=0, ref_mask=None, fit=True, remap=False, mass=False, top=None, dtype='ndarray', frame_indices=None)¶Compute symmetry-corrected RMSD
Parameters: | traj : Trajectory-like mask : str, default ‘’ (all atoms) ref : {int, Frame}, default 0 (first frame) ref_mask : {str, None}, default None
fit : Bool, default True
mass : Bool, default False
remap : Bool, default False
dtype : str, default ‘ndarray’
frame_indices : {None, array-like}, default None
|
---|
Examples
>>> import pytraj as pt
>>> traj = pt.load("TYR.nc", "TYR.parm7")
>>> data = pt.symmrmsd(traj, ref=0)
Notes
-----
versionadded: 1.0.6
pytraj.all_actions.
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.
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
grid_spacing : tuple, grid spacing in X-, Y-, Z-dimensions, require size : {None, tuple}, default None
center : {None, tuple}, default None
buffer : float, default 3.0 Angstrom
centermask : str radscale : float, default 1.36 (to match to VMD calculation)
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.
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.all_actions.
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
min_spacing : {None, float}, default None
byres : bool, default False
|
---|
pytraj.all_actions.
esander
(traj=None, prmtop=None, igb=8, mm_options=None, qm_options=None, dtype='dict', frame_indices=None, top=None)¶energy decomposition by calling libsander
Parameters: | traj : Trajectory-like or iterables that produce Frame
prmtop : str or Structure from ParmEd, default=None, optional
igb : GB model, default=8 (GB-Neck2)
mm_options : InputOptions from sander, default=None, optional
qm_options : InputOptions from sander for QMMM, optional top : pytraj.Topology or str, default=None, optional
dtype : str, {‘dict’, ‘dataset’, ‘ndarray’, ‘dataframe’}, default=’dict’
frame_indices : None or 1D array-like, default None
|
---|---|
Returns: | Dict of energies (to be used with DataFrame) or DatasetList |
Notes
This method does not work with pytraj.pmap when you specify mm_options and qm_options. Use pytraj.pmap_mpi with MPI instead.
Work with pytraj.pmap
:
pt.pmap(pt.esander, traj, igb=8, dtype='dict')
Will NOT work with pytraj.pmap
:
import sander
inp = sander.gas_input(8)
pt.pmap(pt.esander, traj, mm_options=inp, dtype='dict')
Why? Because Python need to pickle each object to send to different cores and Python does not know how to pickle mm_options from sander.gas_input(8).
This works with pytraj.pmap_mpi
because pytraj explicitly create mm_options
in each core without pickling.
Examples
Examples are adapted from $AMBERHOME/test/sanderapi
>>> import pytraj as pt
>>> # GB energy
>>> traj = pt.datafiles.load_ala3()
>>> traj.n_frames
1
>>> data = pt.esander(traj, igb=8)
>>> data['gb']
array([-92.88577683])
>>> data['bond']
array([ 5.59350521])
>>> # PME
>>> import os
>>> from pytraj.testing import amberhome
>>> import sander
>>> topfile = os.path.join(amberhome, "test/4096wat/prmtop")
>>> rstfile = os.path.join(amberhome, "test/4096wat/eq1.x")
>>> traj = pt.iterload(rstfile, topfile)
>>> options = sander.pme_input()
>>> options.cut = 8.0
>>> edict = pt.esander(traj=traj, mm_options=options)
>>> edict['vdw']
array([ 6028.95167558])
>>> # GB + QMMM
>>> topfile = os.path.join(amberhome, "test/qmmm2/lysine_PM3_qmgb2/prmtop")
>>> rstfile = os.path.join(amberhome, "test/qmmm2/lysine_PM3_qmgb2/lysine.crd")
>>> traj = pt.iterload(rstfile, topfile)
>>> options = sander.gas_input(8)
>>> options.cut = 99.0
>>> options.ifqnt = 1
>>> qm_options = sander.qm_input()
>>> qm_options.iqmatoms[:3] = [8, 9, 10]
>>> qm_options.qm_theory = "PM3"
>>> qm_options.qmcharge = 0
>>> qm_options.qmgb = 2
>>> qm_options.adjust_q = 0
>>> edict = pt.esander(traj=traj, mm_options=options, qm_options=qm_options)
>>> edict['bond']
array([ 0.00160733])
>>> edict['scf']
array([-11.92177575])
>>> # passing options to `pytraj.pmap`: need to pass string
>>> from pytraj.testing import get_fn
>>> fn, tn = get_fn('tz2')
>>> traj = pt.iterload(fn, tn)
>>> inp_str = 'mm_options = sander.pme_input()'
>>> edict = pt.pmap(pt.esander, traj, mm_options=inp_str, n_cores=2)
>>> edict['dihedral']
array([ 126.39307126, 127.0460586 , 137.26793522, 125.30521069,
125.25110884, 137.69287326, 125.78280543, 125.14530517,
118.41540102, 128.73535036])
pytraj.all_actions.
lie
(traj, mask, options='', dtype='dict', frame_indices=None, top=None)¶LIE