Source code for StructureBlurrer

## TEMPY is a Python library designed to help the user manipulate and analyse atomic structures and density maps from 3D EM. 
## Copyright (c) 2013 Daven Vasishtan,Irene Farabella, Arun Prasad Pandurangan, Harpal Sahota, Frank Alber and Maya Topf


from EMMap import *
#from ProtRep import *
from ProtRep_Biopy import *

from numpy import array, ndarray, zeros, real
from scipy.fftpack import fftn, ifftn
from scipy.ndimage import fourier_gaussian
from scipy.signal import resample

[docs]class StructureBlurrer: """ A class to generates a density map from a structure instance. """ def __init__(self): pass
[docs] def protMap(self, struct, apix, resolution,filename="None"): """ Returns an Map instance sized and centred based on the atomic structure. Arguments: *apix* Angstroms per pixel for the Map to be outputted. *resolution* Target resolution of the outputted map. *sigma_coeff* Sigma width of the Gaussian used to blur the atomic structure. *filename* output name of the map file. """ # Build empty template map based on the size of the protein and the resolution. extr = struct.get_extreme_values() edge = int(2*resolution/apix)+4 x_size = int((extr[1]-extr[0])/apix)+edge y_size = int((extr[3]-extr[2])/apix)+edge z_size = int((extr[5]-extr[4])/apix)+edge # Origin calculated such that the centre of the map is the centre of mass of the protein. x_origin = struct.CoM.x-(apix*x_size/2.0) y_origin = struct.CoM.y-(apix*y_size/2.0) z_origin = struct.CoM.z-(apix*z_size/2.0) newMap = zeros((z_size, y_size, x_size)) fullMap = Map(newMap, [x_origin, y_origin, z_origin], apix, filename) return fullMap
[docs] def mapGridPosition(self, densMap, atom): """ Returns the index of the nearest pixel to an atom, and atom mass (4 values in list form). Arguments: *densMap* Map instance the atom is to be placed on. *atom* Atom instance. """ origin = densMap.origin apix = densMap.apix box_size = densMap.box_size() x_pos = int(round((atom.x-origin[0])/apix,0)) y_pos = int(round((atom.y-origin[1])/apix,0)) z_pos = int(round((atom.z-origin[2])/apix,0)) #print "grid_pos", x_pos,y_pos,z_pos,atom.x-origin[0], atom.y-origin[1], atom.z-origin[2] #if((box_size[0] > x_pos >= 0) and (box_size[1] > y_pos >= 0) and (box_size[2] > z_pos >= 0)): # return (x_pos, y_pos, z_pos, atom.mass) #MODIFIED BY PAP if((box_size[2] > x_pos >= 0) and (box_size[1] > y_pos >= 0) and (box_size[0] > z_pos >= 0)): return (x_pos, y_pos, z_pos, atom.mass) else: return 0 #this two can be merged and be a unique function that return either the density or 1 #added by PAP
[docs] def make_atom_overlay_map(self, densMap, prot): """ Returns a Map instance with atom masses superposed on it. Arguments: *densMap* an empty (all densities zero) Map instance to superpose the atoms onto. *prot* a Structure instance. """ densMap = densMap.copy() for atom in prot.atomList: pos = self.mapGridPosition(densMap, atom) #print pos if pos: densMap.fullMap[pos[2]][pos[1]][pos[0]] += pos[3] return densMap #ADDED BY PAP
[docs] def make_atom_overlay_map1(self, densMap, prot): """ Returns a Map instance with atom locations recorded on the nearest voxel with a value of 1. Arguments: *densMap* an empty (all densities zero) Map instance to superpose the atoms onto. *prot* a Structure instance. """ densMap = densMap.copy() densMap.fullMap = densMap.fullMap * 0 for atom in prot.atomList: pos = self.mapGridPosition(densMap, atom) #print 'overlay index', pos if pos: densMap.fullMap[pos[2]][pos[1]][pos[0]] = 1 return densMap
[docs] def gaussian_blur(self, prot, resolution, densMap=False, sigma_coeff=0.356, normalise=True,filename="None"): """ Returns a Map instance based on a Gaussian blurring of a protein. Arguments: *prot* the Structure instance to be blurred. *resolution* the resolution, in Angstroms, to blur the protein to. *densMap* False to build a Map with dimensions based on the protein, or a Map instance to be used as a template. *sigma_coeff* the sigma value (multiplied by the resolution) that controls the width of the Gaussian. Default values is 0.356. Other values used : 0.187R corresponding with the Gaussian width of the Fourier transform falling to half the maximum at 1/resolution, as used in Situs (Wriggers et al, 1999); 0.356R corresponding to the Gaussian width at 1/e maximum height equaling the resolution, the default in Chimera (Petterson et al, 2004); 0.425R the fullwidth half maximum being equal to the resolution, as used by FlexEM (Topf et al, 2008); 0.5R the distance between the two inflection points being the same length as the resolution, an option in Chimera (Petterson et al, 2004); 1R where the sigma value simply equal to the resolution, as used by NMFF (Tama et al, 2004). *filename* output name of the map file. """ #densMap= your map if you want to compare prot blurred with an exisiting map. #Daven always use that so that it blurred based on the experiment box if not densMap: densMap = self.protMap(prot, min(resolution/4., 3.5), resolution) print "WARNING: Use StructureBlurrer.gaussian_blur_box() to blured a map with a user defined defined cubic box" #from here till newMap.fullMap*=0 are few line of code that create an empty map with the new A/px of 1 #this replace the make_clash_map(apix) function. they do the job but they need to be replaced with something more rigorous x_s = int(densMap.x_size()*densMap.apix) y_s = int(densMap.y_size()*densMap.apix) z_s = int(densMap.z_size()*densMap.apix) newMap = densMap.resample_by_box_size([z_s, y_s, x_s]) newMap.fullMap *= 0 sigma = sigma_coeff*resolution newMap = self.make_atom_overlay_map(newMap, prot) fou_map = fourier_gaussian(fftn(newMap.fullMap), sigma) newMap.fullMap = real(ifftn(fou_map)) newMap = newMap.resample_by_box_size(densMap.box_size()) if normalise: newMap = newMap.normalise() newMap.filename=filename return newMap #add IF
[docs] def gaussian_blur_box(self, prot, resolution, box_size_x, box_size_y, box_size_z, densMap=False, sigma_coeff=0.356, normalise=True,filename="None"): """ Returns a Map instance based on a Gaussian blurring of a protein. Arguments: *prot* the Structure instance to be blurred. *resolution* the resolution, in Angstroms, to blur the protein to. *box_size_x* x dimension of map box in Angstroms. *box_size_y* y dimension of map box in Angstroms. *box_size_z* z dimension of map box in Angstroms. *densMap* False to build a Map with dimensions based on the protein, or a Map instance to be used as a template. *sigma_coeff* the sigma value (multiplied by the resolution) that controls the width of the Gaussian. Default values is 0.356. Other values used : 0.187R corresponding with the Gaussian width of the Fourier transform falling to half the maximum at 1/resolution, as used in Situs (Wriggers et al, 1999); 0.356R corresponding to the Gaussian width at 1/e maximum height equaling the resolution, the default in Chimera (Petterson et al, 2004); 0.425R the fullwidth half maximum being equal to the resolution, as used by FlexEM (Topf et al, 2008); 0.5R the distance between the two inflection points being the same length as the resolution, an option in Chimera (Petterson et al, 2004); 1R where the sigma value simply equal to the resolution, as used by NMFF (Tama et al, 2004). *filename* output name of the map file. """ if not densMap: densMap = self.protMapBox(prot, 1, resolution, box_size_x, box_size_y, box_size_z, filename) #from here till newMap.fullMap*=0 are few line of code that create an empty map with the new A/px of 1 #this replace the make_clash_map(apix) function. they do the job but they need to be replaced with something more rigorous x_s = int(densMap.x_size()*densMap.apix) y_s = int(densMap.y_size()*densMap.apix) z_s = int(densMap.z_size()*densMap.apix) newMap = densMap.resample_by_box_size([z_s, y_s, x_s]) newMap.fullMap *= 0 sigma = sigma_coeff*resolution newMap = self.make_atom_overlay_map(newMap, prot) fou_map = fourier_gaussian(fftn(newMap.fullMap), sigma) newMap.fullMap = real(ifftn(fou_map)) newMap = newMap.resample_by_box_size(densMap.box_size()) if normalise: newMap = newMap.normalise() return newMap #---BANDPASS FILTERING (NOT WORKING YET)--- add by DV# MAKE them PRIVITA _FUNCT #way of filtering the map using "Fourier-like" but it is too slow so abandon the idea. there are quiker and better way # Bsoft is a better way to go. http://lsbr.niams.nih.gov/bsoft/ # not spend time on it.
def _bandpass_blur(self, atomList, densMap, lopass, lomin, lowid, hipass, hiwid): """ WARNING: BANDPASS FILTERING (NOT WORKING YET) """ pass def _bandpass_mask_gaussian(self, densMap, lopass, lopass_min, lowid, hipass, hiwid): """ WARNING: BANDPASS FILTERING (NOT WORKING YET) """ newMap = densMap.copy()#self.make_empty_map(densMap) centre = (array(newMap.box_size[:])-1)/2.0 from time import time for z in range(newMap.box_size[2]): for y in range(newMap.box_size[1]): for x in range(newMap.box_size[0]): t1 = time() dist = sqrt((x-centre[0])**2 + (y-centre[1])**2 + (z-centre[2])**2) t2 = time() newMap[z][y][x] = self.bandpass_eq_gaussian(dist, lopass, lopass_min, lowid, hipass, hiwid) t3 = time() print t2-t1, t3-t2 return newMap def _bandpass_eq_gaussian(self, dist, lopass, lopass_min, lowid, hipass, hiwid): """ WARNING: BANDPASS FILTERING (NOT WORKING YET) """ lp_max = lopass+lowid hp_min = hipass-hiwid if dist <= lp_max: return lopass_min+(1-lopass_min)*exp(-0.5*((dist-lp_max)/lowid)**2) elif lp_max < dist <= hp_min: return 1.0 else: return exp(-0.5*((dist-hp_min)/hiwid)**2) def _bandpass_test(self, lopass, lopass_min, lowid, hipass, hiwid, l_len): """ WARNING: BANDPASS FILTERING (NOT WORKING YET) """ from time import time start = time() a = zeros([l_len]) for x in range(l_len): a[x] = self.bandpass_eq_gaussian(x, lopass, lopass_min, lowid, hipass, hiwid) end = time() print end-start return a