Cloud#
The cloud submodule includes functions to interact with AV/ACV clouds and interfaces with all other submodules.
- class fretraj.cloud.ACV(grid_acv=None, cv_thickness=0, cv_fraction=0, cloud_xyzqt=None, use_LabelLib=True)[source]#
Reweighted accessible contact volume of a covalently linked fluorophore
- Parameters:
grid_acv (LabelLib.Grid3D, optional=None) – accessible volume with added weight labels for free and contact volume (density label FV: 1.0, density label CV: 2.0); the two volumes are subsequently reweighted
cv_thickness (float, optional=0) – width of the contact volume in Angstrom (default: 0, i.e. no CV is calculated) Note: good first approx.: 2*min(dye_radii)
cv_fraction (float, optional=0) – fraction of dyes that are within the contact volume (e.g. as determined by fluorescence anisotropy)
cloud_xyzq (ndarray, optional=None) – array of x-,y-,z-coordinates and corresponding weights with a shape [n_gridpts(+), 4]
use_LabelLib (bool, optional=True) – make use of LabelLib library to compute FRET values and distances
- _reweight_cv(grid, cv_thickness, cv_fraction)[source]#
Reweight the accessible volume based on contact and free volume
- Parameters:
cv_thickness (float) – width of the contact volume in Angstrom (default: 0, i.e. no CV is calculated)
cv_fraction (float) – relative fraction spent inside the contact volume (i.e. stacked to the surface)
- Returns:
ndarray – one-dimensional array of grid points of length n_gridpts
Notes
A good first approximation of the CV thickness is about two times the smallest dye radius
The CV fraction can be calculated from the residual fluorescence anisotropy
- _tag_volume(grid_1d)[source]#
Assign a tag to the grid values depending on their location in the cloud (1: free volume, 2: contact volume)
- Parameters:
grid_1d (ndarray) – one-dimensional array of grid points length n_gridpts
- Returns:
ndarray – one-dimensional array of length n_gridpts
- static _weight_factor(grid_1d, cv_fraction)[source]#
Calculate the weight of contact volume grid points.
Note
The weight factor corresponds to the factor by which the points of the contact volume (CV, dye trapped on surface) are favored over those belonging the free volume (FV, free dye diffusion). This accounts for the (experimentally) determined fraction of dyes populating the CV.
- Parameters:
grid_1d (ndarray) – one-dimensional array of grid points of length n_gridpts
cv_fraction (float) – fraction of dyes that are within the contact volume
- Returns:
float – weight of each grid point that belongs to the contact volume
- class fretraj.cloud.FRET(volume1, volume2, fret_pair, labels, R_DA=None, verbose=True)[source]#
FRET class of distance and transfer efficiency for a pair of donor and acceptor ACVs
- Parameters:
volume1, volume2 (instances of the Volume class)
fret_pair (str) – Distance key specifying the donor acceptor pair
labels (dict) – dye, linker and setup parameters for the accessible volume calculation
R_DA (ndarray) – donor acceptor distance distribution and associate weights (optional)
verbose (bool)
Examples
>>> ft.Molecule(volume1, volume2)
>>> ft.Molecule(volume1, volume2, use_LabelLib=False)
- classmethod from_volumes(volume_list1, volume_list2, fret_pair, labels, R_DA=None)[source]#
Alternative constructor for the FRET class by reading in a list of donor and acceptor volumes
- Parameters:
volume_list1, volume_list2 (array_like) – lists of Volume instances
fret_pair (str) – identifier for donor acceptor pair
labels (dict) – dye, linker and setup parameters for the accessible volume calculation
R_DA (ndarray) – donor acceptor distance distribution and associate weights (optional)
- save_fret(filename)[source]#
Write the FRET calculation to a json file
- Parameters:
filename (str) – filename for FRET prediction
Examples
>>> obj.save_FRET('parameters.json')
- property values#
Pandas Dataframe of FRET parameters
- Returns:
pandas.DataFrame
- class fretraj.cloud.Trajectory(fret, timestep=None, kappasquare=None)[source]#
Trajectory class of distances and transfer efficiencies from the FRET class
- Parameters:
fret (instance of the FRET class)
timestep (int) – time difference between two frames in picoseconds
- property dataframe#
Pandas Dataframe view of the Trajectory object
- Returns:
pandas.DataFrame – Dataframe with distances and transfer efficiencies of the FRET trajectory
- save_traj(filename, format='csv', units='A', header=True, R_kappa_only=False)[source]#
Save the trajectory as a .csv file
- Parameters:
filename (str) – filename for .csv trajectory
format ({‘csv’, ‘txt’})
units ({‘A’, ‘nm’}, optional=’A’) – distance units (‘A’: Angstroms, ‘nm’: nanometers)
header (bool, default=True) – include header in the output file
R_kappa_only (bool, default=False) – include only time, R_DA and kappasquare columns in the output file (as gmx dyecoupl, Gromacs)
- class fretraj.cloud.Volume(structure, site, labels, cloud_xyzqt=None, verbose=True)[source]#
Class object holding the accessible contact volume of a specific labeling position
- Parameters:
structure (mdtraj.Trajectory) – trajectory of atom coordinates loaded from a pdb, xtc or other file
site (str) – reference identifier for the labeling position
labels (dict) – dye, linker and setup parameters for the accessible volume calculation
cloud_xyzqt (ndarray) – array of x-,y-,z-coordinates with corresponding weights and tags with a shape [n_gridpts(+), 5]
- _dye_acc_surf()[source]#
Calculate dye accessible surface by padding the vdW radius with the thickness of the contact volume
- Returns:
ndarray – array with dimensions (5,n_atoms) of marked coordinates and padded vdW radii (n_atoms = number of atoms in mdtraj.Trajectory)
- static _weighted_quantile_1D(arr_1D, weights, quantile)[source]#
Compute the weighted quantile of a 1D-array
- Parameters:
arr_1D (ndarray) – one-dimensional array
weights (ndarray) – weight of each element of the array
quantile (float) – quantile to compute (median=0.5)
- Returns:
quantile (float)
- property attach_xyz#
Get coordinates of the dye attachment point
- Returns:
ndarray – one-dimensional array of x-,y-,z-coordinates of the attachment point
Examples
>>> avobj.attach_xyz
- calc_acv(use_LabelLib)[source]#
Partition the accessible volume into a free and a contact volume
- Returns:
acv (fretraj.cloud.ACV)
Notes
- The attributes from the LabelLib.Grid3D object are:
discStep : float (the grid spacing)
originXYZ : list (x-/y-/z-coordinate of the grid origin)
shape : list (number of grid points in x-/y-/z-direction)
grid : list (flattened list of grid point values)
- Additional attributes are:
grid_3d : ndarray([shape])
n_gridpts : int
tag : ndarray([1,n_gridpts])
See also
calc_av
Calculate accessible volume
References
Examples
>>> avobj.calc_acv()
- calc_av(use_LabelLib)[source]#
Calculate the dye accessible volume
- Returns:
LabelLib.Grid3D – object containing a 3-dimensional array of grid points and additional attributes (see Notes)
Notes
- The attributes of the av LabelLib.Grid3D object are:
discStep : float (the grid spacing)
originXYZ : list (x-/y-/z-coordinate of the grid origin)
shape : list (number of grid points in x-/y-/z-direction)
grid : list (flattened list of grid point values)
See also
calc_acv
Calculate accessible contact volume
References
Examples
>>> avobj.calc_av()
- classmethod from_attachID(structure, site, labels, attachID)[source]#
Alternative constructor for the Volume class by reading in one or multiple attachment points using the same dye and grid parameters
- Parameters:
structure (mdtraj.Trajectory) – trajectory of atom coordinates loaded from a pdb, xtc or other file
site (str) – reference identifier for the labeling position
labels (dict) – dye, linker and setup parameters for the accessible volume calculation
attachID (int or list) – list of attachment ids on the structure to be used for the ACV calculation
- Returns:
list – list of Volume instances
Examples
>>> struct = md.load_pdb('data/2m23.pdb') >>> labels_global = ft.labeling_params('data/labels_global.json') >>> from_attachID(struct, 0, labelsglobals, [2,929])
- classmethod from_frames(structure, site, labels, frames_mdtraj)[source]#
Alternative constructor for the Volume class by reading in one or multiple frames using the same dye and grid parameters
- Parameters:
structure (mdtraj.Trajectory) – trajectory of atom coordinates loaded from a pdb, xtc or other file
site (str) – reference identifier for the labeling position
labels (dict) – dye, linker and setup parameters for the accessible volume calculation
frames_mdtraj (int or list) – list of frames on the trajectory to be used for the ACV calculation
- Returns:
list – list of Volume instances
- static grid2pts(grid_3d, xyz_min, d_xyz, *args)[source]#
Convert 3D-grid with density values to xyz coordinates with a weight (q)
- Parameters:
grid_3d (ndarray([nx,ny,nz])) – 3-dimensional array of grid points with a shape given by n_xyz
xyz_min (ndarray) – origin coordinates of the grid
d_xyz (ndarray) – grid spacing in x-,y- and z-direction
- Returns:
ndarray – array of x-,y-,z-coordinates with corresponding weights and tags with a shape [n_gridpts(+), 5]
Examples
>>> grid_3d = numpy.zeros((3,3,3)) >>> grid_3d[(1,1,1)] = 1 # FV >>> grid_3d[(1,2,1)] = 10 # CV >>> ft.cloud.Volume.grid2pts(grid_3d, (0,0,0), (1,1,1)) array([[ 0.5, 0.5, 0.5, 1. 1], [ 0.5, 1. , 0.5, 10. 1]])
>>> tag_3d = numpy.zeros((3,3,3)) >>> grid_3d[(1,1,1)] = 1 # FV >>> grid_3d[(1,1,1)] = 2 # CV >>> ft.cloud.Volume.grid2pts(grid_3d, (0,0,0), (1,1,1), tag_3d) array([[ 0.5, 0.5, 0.5, 1. 1], [ 0.5, 1. , 0.5, 10. 2]])
- static mean_pos(cloud_xyzqt)[source]#
Calculate mean dye position of the accessible contact volume
- Parameters:
cloud_xyzq (ndarray) – array of x-,y-,z-coordinates and corresponding weights with a shape [n_gridpts(+), 4]
- Returns:
ndarray – mean position
- property mol_xyzr#
Get coordinates and vdW radii of all selected atoms in the structure
- Returns:
ndarray – array of x-,y-,z-coordinates and VdW radii with a shape [n_atoms, 4]
Examples
>>> avobj.mol_xyzr
- static reshape_grid(grid, shape)[source]#
Convert a 1D-grid to a 3D-Grid
- Parameters:
grid (array_like) – one-dimensional grid array/list of length n_gridpts
shape (list)
- Returns:
ndarray – 3-dimensional array of grid points with a shape given by n_xyz
Examples
>>> grid = numpy.full(27,1) >>> ft.cloud.Volume.reshape_grid(grid, (3,3,3)) array([[[1., 1., 1.], [1., 1., 1.], [1., 1., 1.]], ... [[1., 1., 1.], [1., 1., 1.], [1., 1., 1.]]])
- save_acv(filename, format='pdb', **kwargs)[source]#
Write accessible contact volume to file
- Parameters:
filename (str) – filename for ACV
format ({‘pdb’, ‘xyz’, ‘openDX’}) – file format for ACV
**kwargs (dict) – key value pairs for .xyz files are {write_weights : bool, encode_element : bool}
Notes
write_weights includes the weights of each point as an additional column in the .xyz file
encode_element uses the element column (first col.) to encode the FV (=F) or CV (=C)
Examples
>>> obj.save_acv('cloud.xyz', format='xyz', write_weights=False)
- save_mp(filename, format='plain', units='A')[source]#
Write mean dye position to a file
- Parameters:
filename (str) – filename for the mean position
format ({‘plain’, ‘json’}) – file format (plain .txt or .json)
units ({‘A’, ‘nm’}, optional=’A’) – distance units (‘A’: Angstroms, ‘nm’: nanometers)
Examples
>>> avobj.save_mp('mp.dat')
- fretraj.cloud._create_volume_topology(n_atoms)[source]#
Create a topology for an accessible-contact volume
- Parameters:
n_atoms (int) – number of points in the volume
- fretraj.cloud.acv_subsampling(volume_list, verbose=False)[source]#
Subsample the ACV to obtain identical number of points over all volumes. This allows to create mdtraj.Trajectory objects of the ACV which can be also save as .xtc or .xyz files
- Parameters:
volume_list (array_like) – list of Volume instances
verbose (bool) – be verbose about the number of subsamples
- fretraj.cloud.check_labels(labels, verbose=True)[source]#
Check integrity of parameter dictionary
- Parameters:
labels (dict) – position of labels, dye and grid parameters
verbose (bool, optional=True) – be verbose about missing parameter and their fall-back values
- fretraj.cloud.create_acv_traj(volume_list, verbose=False)[source]#
Create a mdtraj.Trajectory object from a volume list
- Parameters:
volume_list (array_like) – list of Volume instances
verbose (bool) – be verbose about the number of subsamples
- fretraj.cloud.labeling_params(param_file, verbose=True)[source]#
Parse the parameter file to get the configuration settings
- Parameters:
param_file (str) – json-formatted parameter filename
- Returns:
dict – position of labels, dye and grid parameters
- fretraj.cloud.load_acv_traj(filename)[source]#
Load an AV, FV and CV from .xyz or .xtc file
- Parameters:
filename (str) – filename of ACV trajectory
Note
Multi model PDB files cannot be loaded as a mdtraj.Trajectory because they are on purpose not subsampled and thus contain different number of points per volume
- fretraj.cloud.load_obj(filename)[source]#
Load a serialized object from a binary file
- Parameters:
filename (str) – filename of pickle obj
- fretraj.cloud.pipeline_frames(structure, donor_site, acceptor_site, labels, frames, fret_pair)[source]#
Create a pipeline to compute multi-frame donor and acceptor ACVs and calculate a FRET trajectory
- Parameters:
structure (mdtraj.Trajectory) – trajectory of atom coordinates loaded from a pdb, xtc or other file
donor_site (str) – reference identifier for the donor labeling position
acceptor_site (str) – reference identifier for the acceptor labeling position
labels (dict) – dye, linker and setup parameters for the accessible volume calculation
frames_mdtraj (list) – list of frames on the trajectory to be used for the ACV calculation
fret_pair (str) – Distance key specifying the donor acceptor pair
Note
Running a pipeline more memory-efficient than running Volume.from_frames() because the ACVs are stored only until the FRET efficiency is calculated. On the downside, no ACV trajectories(.xtc / .xyz) can be saved.
- fretraj.cloud.printProgressBar(iteration, total, prefix='Progress:', suffix='complete', length=20, fill='█')[source]#
Command line progress bar
- Parameters:
iteration (int) – current iteration
total (int) – total number of iterations
prefix (str, optional=’Progress:’) – string before progress bar
suffix (str, optional=’complete’) – string after progress bar
length (int, optional=20) – length of progress bar
fill (str, optional=’█’) – fill character of progress bar
Examples
>>> printProgressBar(0, n) >>> for i in range(n): doSomething printProgressBar(i + 1, n)
- fretraj.cloud.save_acv_traj(filename, volume_list, format='pdb', verbose=False)[source]#
Save a trajectory of ACVs as a full multi model PDB or a subsampled .xyz or .xtc (with identical number of points over all volumes)
- Parameters:
filename (str) – filename for ACV trajectory
volume_list (array_like) – list of Volume instances
format (str) – file format of the ACV (for ‘pdb’ the full cloud is saved without subsampling; for ‘xyz’ k points in the volumes are randomly picked, where k is the number of points in the smallest volume)
separate_CV (bool) – Save contact and free volume in separate files (choose this if you want to retain the CV information)
- fretraj.cloud.save_labels(filename, labels)[source]#
Write the ACV parameters to a .json file
- Parameters:
filename (str) – filename for label parameters
labels (dict) – position of labels, dye and grid parameters
Examples
>>> obj.save_labels('parameters.json')
- fretraj.cloud.save_mp_traj(filename, volume_list, units='A')[source]#
Save a trajectory of dye mean positions as an xyz file
- Parameters:
filename (str) – filename for mean position trajectory
volume_list (array_like) – list of Volume instances
units ({‘A’, ‘nm’}, optional=’A’) – distance units (‘A’: Angstroms, ‘nm’: nanometers)
- fretraj.cloud.save_obj(filename, obj)[source]#
Save a serialized object to a binary file
- Parameters:
filename (str) – filename for pickle object
obj (serializable object)
- fretraj.cloud.save_structure_traj(filename, structure, frames, format='pdb')[source]#
Save selected frames of a trajectory as multi-model PDB or XTC
- Parameters:
filename (str) – filename for structure trajectory
structure (mdtraj.Trajectory) – trajectory of atom coordinates loaded from a pdb, xtc or other file
format ({‘pdb’, ‘xtc’}) – trajectory file format