## 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
##--Global imports--#
from numpy import ndarray, array, int32, float32, zeros, real, diag, ma, sum as numsum
from random import randrange
from scipy.ndimage.interpolation import rotate, shift, affine_transform, zoom
from scipy.ndimage import laplace
from scipy.ndimage.fourier import fourier_shift
from scipy.fftpack import fftn, ifftn, fftshift
from scipy.signal import resample
from scipy.spatial import KDTree
#from scipy.cluster.vq import kmeans2
import sys, datetime
import struct as binary
#--Local imports--#
from Vector import *
[docs]class Map:
"""
A class representing all information from a density map file.
NOTE: Currently it can only read the mrc format.
"""
def __init__(self, fullMap, origin, apix, filename, header=[]):
"""
Read a map and its parameters in to Map class instance.
*filename*
name of map file.
*origin*
origin co-ordinates of the map (x_origin, y_origin, z_origin).
*apix*
grid spacing of map.
*filename*
filename of the Map instance
NOTE: The *filename* 'build++copy' is reserved for copying of other Map class instances."""
self.header = header
self.origin = origin
self.apix = apix
self.filename = filename
self.fullMap = fullMap
def __repr__(self):
box_size = list(self.box_size())
box_size.reverse()
string1 = 'Obtained from ' + self.filename + '\n'
string2 = 'Origin: '+ str(self.origin) + '\n'
string3 = 'Box size (x,y,z): ' + str(box_size) + '\n'
string4 = 'Grid spacing: ' + str(self.apix) + '\n'
string5 = 'Min, max, mean, std: %.3f, %.3f, %.3f, %.3f \n' %(self.min(), self.max(), self.mean(), self.std())
return string1 + string2 + string3 + string4+string5
[docs] def x_origin(self):
"""
Return the x-coordinate of the origin.
"""
return self.origin[0]
[docs] def y_origin(self):
"""
Return the y-coordinate of the origin.
"""
return self.origin[1]
[docs] def z_origin(self):
"""
Return the z-coordinate of the origin.
"""
return self.origin[2]
[docs] def copy(self):
"""
Return a copy of this Map.
"""
copy = Map(self.fullMap.copy(), self.origin[:], self.apix, self.filename, self.header[:])
return copy
[docs] def getMap(self):
"""
Return the array containing the map data.
"""
return self.fullMap
[docs] def box_size(self):
"""
Return the size of the map array, in ZYX format.
"""
return self.fullMap.shape
[docs] def x_size(self):
"""
Return the x size of the map array.
"""
return self.fullMap.shape[2]
[docs] def y_size(self):
"""
Return the y size of the map array.
"""
return self.fullMap.shape[1]
[docs] def z_size(self):
"""
Return the z size of the map array.
"""
return self.fullMap.shape[0]
[docs] def map_size(self):
"""
Return the size of the array fullMap.
"""
return self.fullMap.size
def __getitem__(self, index):
"""
Allows direct indexing of the map array from the Map instance.
ie. map[2][2][1] instead of map.fullMap[2][2][1]
"""
return self.fullMap[index]
# -- Map modification methods. All of them return a new Map instance. -- #
[docs] def scale_map(self, scaling):
"""
Return a new Map instance scaled by scaling factor.
"""
sc = 1./scaling
c_sh = self.pixel_centre()*(1-sc)
newMap = self.copy()
newMap.fullMap = affine_transform(self.fullMap, diag([sc,sc,sc]), offset=[c_sh.z, c_sh.y, c_sh.x])
#newMap.apix = newMap.apix/scaling
#newMap.origin -= array(newMap.origin)*(self.apix-newMap.apix)/2
#newMap.origin = list(newMap.origin)
return newMap
[docs] def resize_map(self, new_size, centre=False):
"""
Return a Map instance resized.
Arguments:
*new_size*
3-ple (x,y,z) giving the box size.
*centre*
default False
"""
newMap = self.copy()
newMap.fullMap = zeros(new_size)
min_box = [min(x,y) for x,y in zip(newMap.box_size(), self.box_size())]
newMap.fullMap[:min_box[0], :min_box[1], :min_box[2]] = self.fullMap[:min_box[0], :min_box[1], :min_box[2]]
return newMap
#add by IF
[docs] def protMapBox(self, struct, apix, resolution,box_size_x,box_size_y,box_size_z,filename):
"""
Returns a Map instance sized and centered based on the atomic structure.
Arguments:
*struct*
the Structure instance.
*apix*
Angstroms per pixel for the output Map.
*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.
*filename*
output name of the map file.
"""
# Build empty template map based on the size of the protein and the resolution.
x_size = int(box_size_x)
y_size = int(box_size_y)
z_size = int(box_size_z)
# 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 normalise(self):
"""
Return a new Map instance with normalised density values using x_new = x - mean / standard_deviation
"""
newMap = self.copy()
if self.fullMap.std() == 0:
pass
else:
newMap.fullMap = (self.fullMap-self.fullMap.mean())/self.fullMap.std()
return newMap
[docs] def rotate_by_axis_angle(self, x, y, z, angle, CoM, rad=False):
"""
Return a new Map instance rotated around its centre.
Arguments:
*angle*
angle (in radians if rad == True, else in degrees) to rotate map.
*x,y,z*
axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.
*CoM*
centre of mass around which map will be rotated.
"""
m = axis_angle_to_matrix(x, y, z, angle, rad)
# Calculate translation needed to rotate around CoM
newCoM = CoM.matrix_transform(m.T)
offset = CoM-newCoM
# Apply transform
newMap = self.matrix_transform(m, offset.x, offset.y, offset.z)
return newMap
[docs] def rotate_by_euler(self, x, y, z, CoM, rad=False):
"""
Return a new Map instance rotated around pivot given by CoM.
Arguments:
*x,y,z*
Euler angles (in radians if rad == True, else in degrees) used to rotate map.
*CoM*
centre of mass around which map will be rotated.
*x, y, z*
translation in angstroms.
"""
m = euler_to_matrix(x, y, z, rad)
# Calculate translation needed to rotate around CoM
newCoM = CoM.matrix_transform(m.T)
offset = CoM-newCoM
# Apply transform
newMap = self.matrix_transform(m, offset.x, offset.y, offset.z)
return newMap
[docs] def randomise_position(self, max_trans, max_rot, CoM, v_grain=30, rad=False, verbose=False):
###ASK how and who wrote it. WHAT is i (guess structure to generate ensamble)
"""
Randomise the position of this Structure instance.
Arguments:
*max_trans*
translation in angstroms.
*max_rot*
Euler angles (in radians if rad == True, else in degrees) used to rotate map.
*CoM*
centre of mass around which map will be rotated.
*v_grain*
number of point to sample. Default is 30.
"""
t_v = random_vector(-v_grain, v_grain).unit()
r_v = random_vector(-v_grain, v_grain).unit()
j = 1
if max_trans <= 0:
t_dist = 0
else:
t_dist = randrange(1,max_trans)
if max_rot <= 0:
r_ang = 0
else:
r_ang = randrange(1,max_rot)
t_v = t_v.times(t_dist)
self.rotate_by_axis_angle(r_v.x, r_v.y, r_v.z, r_ang, CoM, rad=rad)
if verbose:
print r_v.x, r_v.y, r_v.z, r_ang, t_v.x, t_v.y, t_v.z
return self.translate(t_v.x, t_v.y, t_v.z),r_v.x,r_v.y,r_v.z,t_dist,r_ang
[docs] def rotate_by_matrix(self, mat, CoM, rad=False):
"""
Return a new Map instance rotated around pivot given by CoM.
Arguments:
*mat*
3x3 matrix used to rotate map (in radians if rad == True, else in degrees).
*CoM*
rotation pivot, usually the centre of mass around which map will be rotated.
"""
# Calculate translation needed to rotate around CoM
newCoM = CoM.matrix_transform(mat.T)
offset = CoM-newCoM
# Apply transform
newMap = self.matrix_transform(mat, offset.x, offset.y, offset.z)
return newMap
[docs] def change_origin(self, x_origin, y_origin, z_origin):
"""
Return a new Map instance with changed origin.
Arguments:
*x_origin, y_origin, z_origin*
new co-ordinates of origin.
"""
newMap = self.copy()
newMap.origin = (x_origin, y_origin, z_origin)
return newMap
[docs] def shift_origin(self, x_shift, y_shift, z_shift):
"""
Return a new Map instance with changed origin.
Arguments:
*x_origin, y_origin, z_origin*
new co-ordinates of origin.
"""
newMap = self.copy()
newMap.origin = (self.origin[0]+x_shift, self.origin[1]+y_shift, self.origin[2]+z_shift)
return newMap
[docs] def translate(self, x, y, z):
"""
Return a new Map instance shifted by changing origin.
Arguments:
*x,y,z*
translation in angstroms
"""
sh = array([z/self.apix,y/self.apix,x/self.apix])
newMap = self.copy()
newMap.fullMap = shift(newMap.fullMap, sh, cval=self.min())
#f_map = fftn(newMap.fullMap)
#newMap.fullMap = real(ifftn((fourier_shift(f_map, shift))))
return newMap
[docs] def origin_change_maps(self,MapRef):
"""
Return a new Map instance with origin changed accordingly to Reference Map (*MapRef*)
"""
newMap = self.copy()
origin_shift = [y-x for x,y in zip(newMap.origin, MapRef.origin)]
#m2 = m2.shift_origin(origin_shift[0],origin_shift[1],origin_shift[2])
newMap.translate(origin_shift[0],origin_shift[1],origin_shift[2])
newMap.origin = MapRef.origin[:]
return newMap
[docs] def threshold_map(self, minDens, maxDens):
"""
Return a new Map instance containing only density values between the specified min and max.
Arguments:
*minDens*
minimum density threshold
*maxDens*
maximum density threshold
"""
a = self.fullMap.copy()
a = a*(a<maxDens)*(a>minDens)
b = self.copy()
b.fullMap = a
return b
#add by AP
#for SS_CCF score
# maybe change their name
def _get_maskArray(self, densValue):
"""ADDED by APP to use with SSC_CCC score"""
mask_array = ma.masked_less_equal(self.fullMap, densValue)
return ma.getmask(mask_array)
def _get_maskMap(self, maskArray):
"""ADDED by APP to use with SSC_CCC score"""
return ma.masked_array(self.fullMap, mask=maskArray)
[docs] def make_bin_map(self, cutoff):
"""
Return a new Map instance that has been binarised. All densities above the cutoff take a value of '1', else '0'.
Arguments:
*cutoff*
cutoff density value
"""
binMap = self.copy()
newMap = self.fullMap > cutoff
binMap.fullMap = newMap*-1
return binMap
[docs] def make_clash_map(self,apix=1.0):
#note this is fine for the structure blured
# look at note there DAVE NEED TO CTRL
## that should be used in the blured but must be checked
"""
NOTE: NEEED TO BE CHECKED.
Return an empty Map instance with set Angstrom per pixel sampling (default is 1)
"""
grid = zeros((self.box_size()[2]*self.apix/apix, self.box_size()[1]*self.apix/apix, self.box_size()[0]*self.apix/apix))
clashMap = self.copy()
clashMap.fullMap = grid
clashMap.apix = apix
return clashMap
[docs] def resample_by_apix(self, new_apix):
"""
Return a Map instance instance resampled based on new_apix.
"""
new_map = self.copy()
apix_ratio = self.apix/new_apix
#new_map.apix = new_apix
new_map.fullMap = resample(new_map.fullMap, self.z_size()*apix_ratio, axis=0)
new_map.fullMap = resample(new_map.fullMap, self.y_size()*apix_ratio, axis=1)
new_map.fullMap = resample(new_map.fullMap, self.x_size()*apix_ratio, axis=2)
new_map.apix = (self.apix*self.box_size()[2])/new_map.box_size()[2]
return new_map
[docs] def resample_by_box_size(self, new_box_size):
"""
Return a Map instance instance resampled based on new_box_size.
"""
new_map = self.copy()
new_map.fullMap = resample(new_map.fullMap, new_box_size[0], axis=0)
new_map.fullMap = resample(new_map.fullMap, new_box_size[1], axis=1)
new_map.fullMap = resample(new_map.fullMap, new_box_size[2], axis=2)
new_map.apix = (self.apix*self.box_size()[2])/new_box_size[2]
return new_map
# -- Map modifications involving filtering. All still return a new Map instance -- #
[docs] def fourier_filter(self, resolution):
"""
NOTE: UNDER CONSTRUCTION
Return a Map instance of the density map filtered to specific resolution.
Arguments:
*resolution*
highest frequency value to allow into filtered map
"""
a = self.fourier_transform()
x = self.apix*self.box_size()[0]/(2*resolution)
y = self.apix*self.box_size()[1]/(2*resolution)
z = self.apix*self.box_size()[2]/(2*resolution)
x0 = a.x_size()/2.0 + 0.5
y0 = a.y_size()/2.0 + 0.5
z0 = a.z_size()/2.0 + 0.5
#print x,y,z
#print sum(a.fullMap)
for k in range(a.z_size()/2):
for j in range(a.y_size()/2):
for i in range(a.x_size()/2):
if (k/z)**2 + (j/y)**2 + (i/x)**2 >= 1:
a.fullMap[k][j][i] = 0
a.fullMap[k][j+a.y_size()/2][i] = 0
a.fullMap[k][j][i+a.x_size()/2] = 0
a.fullMap[k][j+a.y_size()/2][i+a.x_size()/2] = 0
a.fullMap[k+a.z_size()/2][j][i] = 0
a.fullMap[k+a.z_size()/2][j+a.y_size()/2][i] = 0
a.fullMap[k+a.z_size()/2][j][i+a.x_size()/2] = 0
a.fullMap[k+a.z_size()/2][j+a.y_size()/2][i+a.x_size()/2] = 0
#print sum(a.fullMap)
#===============================================================================
# """a.fullMap[z:a.z_size/2] *= 0
# a.fullMap[(a.z_size/2):(a.z_size/2)+z] *= 0
#
# a.fullMap[:,y:a.y_size/2] *= 0
# a.fullMap[:,(a.y_size/2):(a.y_size/2)+y] *= 0
#
# a.fullMap[:,:,x:a.x_size/2] *= 0
# a.fullMap[:,:,(a.x_size/2):(a.x_size/2)+x] *= 0"""
#===============================================================================
b = ifftn(a.fullMap)
c = self.copy()
c.fullMap = b.copy()
return c
[docs] def laplace_filtered(self):
"""
Return a Map instance of a Laplace filtered density map.
"""
m = self.copy()
m.fullMap = laplace(self.fullMap)
return m
# -- Methods that obtain information from the density map -- #
[docs] def get_vectors(self):
"""
Return an array of all non-zero density points in the form of Vector instances.
"""
a = []
for z in range(len(self.fullMap)):
for y in range(len(self.fullMap[z])):
for x in range(len(self.fullMap[z][y])):
if self.fullMap[z][y][x] != 0:
a.append((Vector((x*self.apix)+self.origin[0], (y*self.apix)+self.origin[1], (z*self.apix)+self.origin[2]),self.fullMap[z][y][x]))
return array(a)
[docs] def get_line_between_points(self, p1, p2):
"""
Return an array of float values representing a line of density values between two points on the map.
Arguments:
*p1, p2*
Vector instances of the end points co-ordinates of the line."""
v1 = p1.minus(Vector(self.origin[0], self.origin[1], self.origin[2])).times(1.0/self.apix)
v2 = p2.minus(Vector(self.origin[0], self.origin[1], self.origin[2])).times(1.0/self.apix)
v3 = v2.minus(v1)
noOfPoints = int(round(2*v3.mod()/self.apix))
points = []
for x in range(0, noOfPoints+1):
p = v1.plus(v3.times(float(x)/noOfPoints))
points.append(self.fullMap[p.z][p.y][p.x])
return array(points)
[docs] def get_com(self):
"""
Return Vector instance of the centre of mass of the map.
"""
total_x_moment = 0.0
total_y_moment = 0.0
total_z_moment = 0.0
total_mass = 0.0
min_mass = self.min()
vectorMap = self.get_vectors()
for v in vectorMap:
m = v[1]+min_mass
total_mass += m
total_x_moment += m*v[0].x
total_y_moment += m*v[0].y
total_z_moment += m*v[0].z
x_CoM = total_x_moment/total_mass
y_CoM = total_y_moment/total_mass
z_CoM = total_z_moment/total_mass
return Vector(x_CoM, y_CoM, z_CoM)
[docs] def pixel_centre(self):
"""
Return Vector instance of the centre of the map in pixels.
"""
x_centre = self.x_size()/2
y_centre = self.y_size()/2
z_centre = self.z_size()/2
return Vector(x_centre, y_centre, z_centre)
[docs] def centre(self):
"""
Return Vector instance of the centre of the map in Angstroms.
"""
x_centre = self.origin[0]+(self.apix*self.x_size()/2)
y_centre = self.origin[1]+(self.apix*self.y_size()/2)
z_centre = self.origin[2]+(self.apix*self.z_size()/2)
return Vector(x_centre, y_centre, z_centre)
[docs] def mean(self):
"""
Return mean density value of map.
"""
return self.fullMap.mean()
[docs] def std(self):
"""
Return standard deviation of density values in map.
"""
return self.fullMap.std()
[docs] def min(self):
"""
Return minimum density value of the map.
"""
return self.fullMap.min()
[docs] def max(self):
"""
Return maximum density value of the map.
"""
return self.fullMap.max()
[docs] def vectorise_point(self, x, y, z):
"""
Return a tuple of the Angstrom co-ordinates and density value of a particular density point in map.
Arguments:
*x, y, z*
co-ordinates of the density point to be vectorised.
"""
v_x = (self.apix*x)+self.origin[0]
v_y = (self.apix*y)+self.origin[1]
v_z = (self.apix*z)+self.origin[2]
dens = self.fullMap[z][y][x]
return Vector(v_x, v_y, v_z)
[docs] def get_significant_points(self):
"""
Return an array of tuples of the map co-ordinates and density values of all points with a density
greater than the mean + standard deviation.
"""
sig_points = []
boo = self.fullMap > (self.fullMap.mean() + self.fullMap.std())
for z in range(self.z_size()):
for y in range(self.y_size()):
for x in range(self.x_size()):
if boo[z][y][x]:
sig_points.append(array([z,y,x,self.fullMap[z][y][x]]))
return array(sig_points)
def _get_random_significant_pairs(self, amount):
"""
Return an array of tuple pairs of significant points randomly chosen from 'get_significant_points' function.
For use with the DCCF and DLSF scoring functions.
Arguments:
*amount*
number of significant point pairs to return.
"""
sig_points = self.get_significant_points()
sig_pairs = []
size = len(sig_points)
for r in range(amount):
fst = sig_points[randrange(size)]
snd = sig_points[randrange(size)]
new_value = array([fst[0], fst[1], fst[2], snd[0], snd[1], snd[2], fst[3]-snd[3]])
sig_pairs.append(new_value)
return array(sig_pairs)
[docs] def makeKDTree(self, minDens, maxDens):
"""
Returns k-dimensional tree of points in the map with values between those chosen.
Arguments:
*minDens*
minimum density value to include in k-dimensional tree.
*maxDens*
maximum density value to include in k-dimensional tree.
"""
maplist = self.get_pos(minDens, maxDens)
return KDTree(maplist)
[docs] def get_pos(self, minDens, maxDens):
"""
Returns array of 3-ples of points in the map with values between those chosen.
Arguments:
*minDens*
minimum density value to include in array.
*maxDens*
maximum density value to include in array.
"""
a = []
for z in range(len(self.fullMap)):
for y in range(len(self.fullMap[z])):
for x in range(len(self.fullMap[z][y])):
if (self.fullMap[z][y][x] > minDens) and (self.fullMap[z][y][x] < maxDens):
a.append((x,y,z))
return array(a)
[docs] def get_normal_vector(self, x_pos, y_pos, z_pos):
"""
Return normal vector at point specified.
Point calculated using 3SOM algorithm used by Ceulemans H. & Russell R.B. (2004).
Arguments:
*x_pos, y_pos, z_pos*
pixel in map on which to calculate normal vector.
"""
internal_vecs = []
for x in range(x_pos-1, x_pos+2):
for y in range(y_pos-1, y_pos+2):
for z in range(z_pos-1, z_pos+2):
if (x_pos, y_pos, z_pos) == (x,y,z):
pass
else:
if self.fullMap[z][y][x] > self.fullMap[z_pos][y_pos][x_pos]:
internal_vecs.append(Vector(x-x_pos,y-y_pos,z-z_pos).unit())
sub_vector = Vector(0,0,0)
for v in internal_vecs:
sub_vector = sub_vector+v
return sub_vector.unit()
[docs] def represent_normal_vectors(self, min_threshold, max_threshold):
"""
Return Structure instance representing normal vectors of density points specified.
Arguments:
*min_threshold, max_threshold*
minimum/maximum values to include in normal vector representation.
"""
atomList = []
template = 'HETATM 1 C NOR A 1 23.161 39.732 -25.038 1.00 10.00 C'
m = self.copy()
print m.origin
for x in range(1, (m.box_size()[0]-1)):
for y in range(1, (m.box_size()[1]-1)):
for z in range(1, (m.box_size()[2]-1)):
if m.fullMap[z][y][x] > min_threshold and m.fullMap[z][y][x] < max_threshold:
#n_vec = m.get_normal_vector(x,y,z, min_threshold)
n_vec = m.get_normal_vector(x,y,z)
n_vec = n_vec.unit()
pos_vec = Vector((x*m.apix)+m.origin[0], (y*m.apix)+m.origin[1], (z*m.apix)+m.origin[2])
#a = BioPyAtom(template)
a = template.BioPyAtom()
a.x = pos_vec.x
a.y = pos_vec.y
a.z = pos_vec.z
b = BioPyAtom(template)
b.x = pos_vec.x + n_vec.x
b.y = pos_vec.y + n_vec.y
b.z = pos_vec.z + n_vec.z
c = BioPyAtom(template)
c.x = pos_vec.x + 2*n_vec.x
c.y = pos_vec.y + 2*n_vec.y
c.z = pos_vec.z + 2*n_vec.z
atomList.append(a)
atomList.append(b)
atomList.append(c)
try:
s = BioPy_Structure(atomList)
except ZeroDivisionError:
return atomList
s.renumber_atoms()
#for x in range(1, len(atomList),2):
# s.footer += 'CONECT '+str(x)+' '+str(x+1)+'\n'
return s
def _set_min_zero(self):
newMap = self.copy()
if newMap.min()<0:
print "here"
low=(0-newMap.min())
print low
newMap.fullMap=newMap.fullMap+abs(low)
return newMap
else:
return newMap
#===============================================================================
#
# #These two values should be calculated for the experimental map, and only
# need to be calculated once, at the beginning. Also note that when
# running the normal_vector_score, the experimental map should be the map1
# parameter (ie, the first one)
#===============================================================================
[docs] def get_min_threshold(self, molWeight, low, high):
"""
Return minimum density value.
Volume of pixels with greater density than output is equivalent to volume given by molecular weight of protein.
Uses recursive algorithm.
Arguments:
*molWeight*
molecular weight of protein, using get_prot_mass_from_atoms().
*low, high*
minimum and maximum values between which the threshold will be taken.
Initial values should be given by minimum and maximum density values in map.
"""
#print low,high
if high-low < 0.0000002 or high < low:
est_weight = long(numsum(self.fullMap > low)*self.apix**3/1210)
print 'Exact molecular weight cannot be found. Approx. weight of '+str(est_weight)+' used instead.'
return low
thr = low+(high-low)/2
this_weight = long(numsum(self.fullMap > thr)*self.apix**3/1210)
#this_weight = long((self.fullMap > thr).sum()*self.apix**3/1210)
#print thr, this_weight, this_weight.sum()
if this_weight == long(molWeight):
return thr
elif this_weight > molWeight:
return self.get_min_threshold(molWeight, thr, high)
elif this_weight < molWeight:
return self.get_min_threshold(molWeight, low, thr)
[docs] def get_max_threshold(self, min_thr, noOfPoints, low, high, err_percent=1):
"""
Return maximum density value.
Given minimum threshold value, will calculate upper bound density value which gives specified number of points with density values between the two.
Uses recursive algorithm.
Arguments:
*min_thr*
minimum threshold, normally given by get_min_threshold method.
*noOfPoints*
Number of points to use in the normal vector score - try first with 10% (maybe 5%better) of the number of points in the map ( round((self.map_size())*0.1)
*low, high*
minimum and maximum values between which the threshold will be taken.
Initial values should be given by *min_thr parameter* and maximum density values in map.
*err_percent*
default value of 1. Allowed to find a max threshold that include a 1% error.
"""
#print low, high
if high-low < 0.0000002 or high<low:
est_weight = numsum((self.fullMap < low)*(self.fullMap > min_thr))
print 'Cannot find value to match number of pixels. Try increasing size of err_percent'
return 0
thr = low+(high-low)/2
this_no_of_points = numsum((self.fullMap < thr)*(self.fullMap > min_thr))
#print thr, this_no_of_points
if abs(this_no_of_points - noOfPoints) < err_percent*noOfPoints/100.:
return thr
elif this_no_of_points < noOfPoints:
return self.get_max_threshold(min_thr, noOfPoints, thr, high)
elif this_no_of_points > noOfPoints:
return self.get_max_threshold(min_thr, noOfPoints, low, thr)
# DO THIS (SHRINK BOX SIZE TO AROUND SD VALUE)
def _shrink_map(self, sd=2.):
pass
# -- Writing Map instance to files -- #
def _write_to_xplor_file(self, xplorFileName):
"""OBSOLETE
xplorFileName = name of file to write to. Note that this function does not automatically append a .xplor suffix."""
xplor = '\n 1 !NTITLE\n'
xplor += 'REMARKS '+'"' + xplorFileName + '"' + ' written by ME!\n'
xplor += self._pad_grid_line_no(self.z_size()) + self._pad_grid_line_no(0) + self._pad_grid_line_no(self.z_size()-1) + \
self._pad_grid_line_no(self.y_size()) + self._pad_grid_line_no(0) + self._pad_grid_line_no(self.y_size()-1) + \
self._pad_grid_line_no(self.x_size()) + self._pad_grid_line_no(0) + self._pad_grid_line_no(self.x_size()-1) + '\n'
xplor += self._convert_point_to_string(self.apix*self.z_size()) + self._convert_point_to_string(self.apix*self.y_size()) + \
self._convert_point_to_string(self.apix*self.x_size())
xplor += self._convert_point_to_string(90.0) + self._convert_point_to_string(90.0) + self._convert_point_to_string(90.0) + '\n'
xplor += 'ZYX\n'
flatList = self.fullMap.flatten()
blockSize = self.z_size()*self.y_size()
blockNo = 0
offset = 0
for point in range(len(flatList)):
if ((point-offset)%6 == 0 and point%blockSize != 0):
xplor += '\n'
if point%blockSize == 0:
xplor += '\n'+ self._pad_grid_line_no(blockNo) +'\n'
blockNo += 1
offset = point%6
xplor += self._convert_point_to_string(real(flatList[point]))
xplor += '\n -9999'
f = file(xplorFileName, 'w')
f.write(xplor)
f.close()
def _write_to_situs_file(self, situsFileName):
"""One day I will do this."""
pass
[docs] def write_to_MRC_file(self, mrcfilename):
"""
writing to MRC file
"""
h = self.header_to_binary()
maparray = array(self.fullMap, dtype='float32')
f = open(mrcfilename, 'wb')
f.write(h)
f.write(maparray.tostring())
f.close()
def _pad_grid_line_no(self, no):
"""Private function to help write data to map files."""
s = str(no)
spaces = ''
for x in range(8-len(s)):
spaces += ' '
s = spaces + s
return s
def _convert_point_to_string(self, point):
"""Private function to help write data to map files."""
exp = 0
sign = ''
if(abs(point)<0.0001):
point = 0.0
if point >= 0:
sign = '+'
else:
sign = '-'
while abs(point) >= 10.0:
point /= 10.0
exp += 1
pointString = str(point)
if len(pointString) < 7:
for x in range(len(pointString), 7):
pointString += '0'
elif len(pointString) > 7:
pointString = pointString[:7]
pointString += 'E' + sign + '0' + str(exp)
return ' '+pointString
def _get_component_volumes(self,struct, apix, blurrer):
"""Private function for check on Map instance."""
mapCoM = self.get_com()
ssplit = struct.split_into_chains()
temp_grid = self.make_clash_map(apix)
overlay_maplist = []
cvol = []
for x in ssplit:
tx = mapCoM[0] - x.CoM[0]
ty = mapCoM[1] - x.CoM[1]
tz = mapCoM[2] - x.CoM[2]
x.translate(tx,ty,tz)
overlay_maplist.append(blurrer.make_atom_overlay_map1(temp_grid, x))
for y in overlay_maplist:
cvol.append(y.fullMap.sum()*(apix**3))
return cvol
#added by PAP
[docs] def map_rotate_by_axis_angle(self, x, y, z, angle, CoM, rad=False):
"""
Return a new Map instance rotated around its centre.
Arguments:
*angle*
angle (in radians if rad == True, else in degrees) to rotate map.
*x,y,z*
axis to rotate about, ie. x,y,z = 0,0,1 rotates the Map round the xy-plane.
*CoM*
centre of mass around which map will be rotated.
"""
m = axis_angle_to_matrix(x, y, z, angle, rad)
#print Vector(x,y,z).unit()
# Calculate translation needed to rotate around CoM
newCoM = CoM.matrix_transform(m.T)
offset = CoM-newCoM
# Apply transform
newMap = self.matrix_transform(m, offset.x, offset.y, offset.z)
return newMap
#added by PAP
[docs] def map_randomise_position(self, i, max_trans, max_rot, CoM, v_grain=30, rad=False, verbose=False):
"""
Randomise the position of the map instance and return the modified map.
Arguments:
*ARUN need to add it*
"""
t_v = random_vector(-v_grain, v_grain).unit()
r_v = random_vector(-v_grain, v_grain).unit()
j = 1
if max_trans <= 0:
t_dist = 0
else:
t_dist = randrange(1,max_trans)
if max_rot <= 0:
r_ang = 0
else:
r_ang = randrange(1,max_rot)
t_v = t_v.times(t_dist)
#self.rotate_by_axis_angle(r_v.x, r_v.y, r_v.z, r_ang, t_v.x, t_v.y, t_v.z, rad=rad)
rotmap = self.rotate_by_axis_angle(r_v.x, r_v.y, r_v.z, r_ang, CoM, rad=rad)
if verbose:
print '%6d %7.3f %7.3f %7.3f %5.2f %5.2f' % (i,r_v.x,r_v.y,r_v.z,t_dist,r_ang)
return rotmap.translate(t_v.x, t_v.y, t_v.z),r_v.x,r_v.y,r_v.z,t_dist,r_ang