Package mmLib :: Module Superposition
[hide private]
[frames] | no frames]

Source Code for Module mmLib.Superposition

  1  ## Copyright 2002-2010 by PyMMLib Development Group (see AUTHORS file) 
  2  ## This code is part of the PyMMLib distribution and governed by 
  3  ## its license.  Please see the LICENSE file that should have been 
  4  ## included as part of this package. 
  5  """Class for least-squares structural superposition.  Uses a quaternion 
  6  method which avoids improper rotations. 
  7  """ 
  8  import math 
  9   
 10  try: 
 11      import numpy 
 12      try: 
 13          from numpy.oldnumeric import linear_algebra as linalg 
 14      except ImportError: 
 15          from numpy.linalg import old as linalg 
 16  except ImportError: 
 17      import NumericCompat as numpy 
 18      from NumericCompat import linalg 
 19       
 20  import AtomMath 
 21   
22 -def QuaternionToRotationMatrix(q):
23 """Create a rotation matrix from q quaternion rotation. 24 Quaternions are typed as Numeric Python arrays of length 4. 25 """ 26 assert numpy.allclose(math.sqrt(numpy.dot(q,q)), 1.0) 27 28 q0, q1, q2, q3 = q 29 30 b0 = 2.0*q0 31 b1 = 2.0*q1 32 b2 = 2.0*q2 33 b3 = 2.0*q3 34 35 q00 = b0*q0-1.0 36 q01 = b0*q1 37 q02 = b0*q2 38 q03 = b0*q3 39 40 q11 = b1*q1 41 q12 = b1*q2 42 q13 = b1*q3 43 44 q22 = b2*q2 45 q23 = b2*q3 46 47 q33 = b3*q3 48 49 return numpy.array([ [q00+q11, q12-q03, q13+q02], 50 [q12+q03, q00+q22, q23-q01], 51 [q13-q02, q23+q01, q00+q33] ], float)
52 53
54 -class SuperpositionResults(object):
55 """Returns the results of a superposition. 56 """
57 - def __init__(self, quaternion, source_origin, destination_origin, rmsd, num_atoms):
58 self.Q = quaternion 59 self.R = QuaternionToRotationMatrix(quaternion) 60 self.src_origin = source_origin 61 self.dst_origin = destination_origin 62 self.rmsd = rmsd 63 self.num_atoms = num_atoms
64
65 - def transform(self, position):
66 """Transforms a source position to its aligned position. 67 """ 68 position = position - self.src_origin 69 position = numpy.dot(self.R, position) 70 position = position + self.dst_origin 71 return position
72 73
74 -def SuperimposePoints(src_points, dst_points):
75 """Takes two 1:1 set of points and returns a 3x3 rotation matrix and 76 translation vector. 77 """ 78 num_points = src_points.shape[0] 79 80 ## shift both sets of coordinates to their centroids 81 src_org = numpy.add.reduce(src_points) / float(src_points.shape[0]) 82 dst_org = numpy.add.reduce(dst_points) / float(dst_points.shape[0]) 83 84 X = numpy.add(src_points, -src_org) 85 Y = numpy.add(dst_points, -dst_org) 86 87 xy2n = 0.0 88 89 R = numpy.zeros((3,3), float) 90 91 for k in xrange(num_points): 92 x = X[k] 93 y = Y[k] 94 95 xy2n += numpy.add.reduce(x*x) + numpy.add.reduce(y*y) 96 97 R[0,0] += x[0]*y[0] 98 R[0,1] += x[0]*y[1] 99 R[0,2] += x[0]*y[2] 100 101 R[1,0] += x[1]*y[0] 102 R[1,1] += x[1]*y[1] 103 R[1,2] += x[1]*y[2] 104 105 R[2,0] += x[2]*y[0] 106 R[2,1] += x[2]*y[1] 107 R[2,2] += x[2]*y[2] 108 109 F = numpy.zeros((4,4), float) 110 F[0,0] = R[0,0] + R[1,1] + R[2,2] 111 F[0,1] = R[1,2] - R[2,1] 112 F[0,2] = R[2,0] - R[0,2] 113 F[0,3] = R[0,1] - R[1,0] 114 115 F[1,0] = F[0,1] 116 F[1,1] = R[0,0] - R[1,1] - R[2,2] 117 F[1,2] = R[0,1] + R[1,0] 118 F[1,3] = R[0,2] + R[2,0] 119 120 F[2,0] = F[0,2] 121 F[2,1] = F[1,2] 122 F[2,2] =-R[0,0] + R[1,1] - R[2,2] 123 F[2,3] = R[1,2] + R[2,1] 124 125 F[3,0] = F[0,3] 126 F[3,1] = F[1,3] 127 F[3,2] = F[2,3] 128 F[3,3] =-R[0,0] - R[1,1] + R[2,2] 129 130 evals, evecs = linalg.eigenvectors(F) 131 132 i = numpy.argmax(evals) 133 eval = evals[i] 134 evec = evecs[i] 135 136 msd = (xy2n - 2.0*eval) / num_points 137 if msd < 0.0: 138 rmsd = 0.0 139 else: 140 rmsd = math.sqrt(msd) 141 142 return SuperpositionResults(evec, src_org, dst_org, rmsd, num_points)
143 144
145 -def SuperimposePositions(position_tuple_list):
146 """Superimposes a list of 2-tuple atom pairs. 147 """ 148 n = len(position_tuple_list) 149 a1 = numpy.zeros((n,3), float) 150 a2 = numpy.zeros((n,3), float) 151 152 for i in xrange(n): 153 pos1, pos2 = position_tuple_list[i] 154 155 a1[i] = pos1 156 a2[i] = pos2 157 158 return SuperimposePoints(a1, a2)
159 160
161 -def SuperimposeAtoms(atom_pair_list):
162 """Superimposes a list of 2-tuple atom pairs. 163 """ 164 n = len(atom_pair_list) 165 a1 = numpy.zeros((n,3), float) 166 a2 = numpy.zeros((n,3), float) 167 168 for i in xrange(n): 169 atm1, atm2 = atom_pair_list[i] 170 171 a1[i] = atm1.position 172 a2[i] = atm2.position 173 174 return SuperimposePoints(a1, a2)
175 176
177 -def SuperimposeAtomsOutlierRejection(alist, rmsd_cutoff = 1.0, max_cycles = 100):
178 """Superimpose two homologous protein chains. The argument alist is a list of 179 2-tuples. The 2-tuples are the 1:1 atoms to superimpose. The alignment 180 procedure incrementally omits atoms with large deviations until the rmsd of 181 the least squares superposition is less than or equal to rmsd_cutoff, or the 182 number of cycles exceeds max_cycles. 183 """ 184 for cycle in xrange(max_cycles): 185 sresult = SuperimposeAtoms(alist) 186 print "Cycle %d NA: %5d RMSD %6.2f" % (cycle, sresult.num_atoms, sresult.rmsd) 187 if sresult.rmsd <= rmsd_cutoff: 188 return sresult 189 190 deviation = [] 191 for i in xrange(len(alist)): 192 satm, datm = alist[i] 193 spos = sresult.transform(satm.position) 194 deviation.append((AtomMath.length(spos - datm.position), i)) 195 deviation.sort() 196 197 outliers = deviation[-10:] 198 outliersr = [(x[1], x[0]) for x in outliers] 199 outliersr.sort() 200 outliersr.reverse() 201 202 for i, d in outliersr: 203 if d < rmsd_cutoff: continue 204 del alist[i] 205 206 return None
207 208 ## <testing>
209 -def test_module():
210 import random 211 import AtomMath 212 import FileIO 213 import Structure 214 215 R = AtomMath.rmatrixu(numpy.array((0.0, 0.0, 1.0), float), math.pi/2.0) 216 217 struct1 = FileIO.LoadStructure(fil="8rxn.pdb") 218 struct2 = FileIO.LoadStructure(fil="8rxn.pdb") 219 220 chn1 = struct1.get_chain("A") 221 chn2 = struct2.get_chain("A") 222 223 rc = lambda: 0.1 * (random.random() - 1.0) 224 for atm in chn2.iter_atoms(): 225 atm.position = numpy.dot(R, atm.position) + numpy.array((rc(),rc(),rc()),float) 226 227 alist = [] 228 229 for atm1 in chn1.iter_atoms(): 230 if atm1.name != "CA": 231 continue 232 233 atm2 = chn2.get_equivalent_atom(atm1) 234 if atm2 == None: continue 235 236 alist.append((atm1, atm2)) 237 238 sup = SuperimposeAtoms(alist) 239 240 R = sup.R 241 Q = sup.Q 242 print Q 243 print R 244 245 so = sup.src_origin 246 do = sup.dst_origin 247 248 sup1 = Structure.Structure(structure_id = "JMP1") 249 for atm in chn1.iter_atoms(): 250 atm.position = numpy.dot(R, atm.position - so) 251 sup1.add_atom(atm) 252 FileIO.SaveStructure(fil="super1.pdb", struct=sup1) 253 254 sup2 = Structure.Structure(structure_id = "JMP2") 255 for atm in chn2.iter_atoms(): 256 atm.position = atm.position - do 257 sup2.add_atom(atm) 258 FileIO.SaveStructure(fil="super2.pdb", struct=sup2)
259 260 if __name__ == "__main__": 261 test_module() 262 ## </testing> 263