00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019 """
00020 Wrappers for the tmalign and tmscore utilities.
00021
00022 References:
00023
00024 tmscore: Yang Zhang and Jeffrey Skolnick, Proteins 2004 57: 702-710
00025 tmalign: Y. Zhang and J. Skolnick, Nucl. Acids Res. 2005 33, 2302-9
00026
00027
00028 Authors: Pascal Benkert, Marco Biasini
00029 """
00030
00031 import subprocess, os, tempfile, platform
00032 from ost import settings, io, geom
00033
00034 def _SetupFiles(models):
00035
00036 tmp_dir_name=tempfile.mkdtemp()
00037 for index, model in enumerate(models):
00038 io.SavePDB(model, os.path.join(tmp_dir_name, 'model%02d.pdb' % (index+1)))
00039 return tmp_dir_name
00040
00041 def _CleanupFiles(dir_name):
00042 import shutil
00043 shutil.rmtree(dir_name)
00044
00045 class TMAlignResult:
00046 def __init__(self, rmsd, tm_score, aligned_length, transform):
00047 self.rmsd=rmsd
00048 self.tm_score=tm_score
00049 self.aligned_length=aligned_length
00050 self.transform=transform
00051
00052 def _ParseTmAlign(lines):
00053 info_line=lines[11].split(',')
00054 aln_length=float(info_line[0].split('=')[1].strip())
00055 rmsd=float(info_line[1].split('=')[1].strip())
00056 tm_score=float(info_line[2].split('=')[1].strip())
00057 tf1=[float(i.strip()) for i in lines[15].split()]
00058 tf2=[float(i.strip()) for i in lines[16].split()]
00059 tf3=[float(i.strip()) for i in lines[17].split()]
00060 rot=geom.Mat3(tf1[2], tf1[3], tf1[4], tf2[2], tf2[3],
00061 tf2[4], tf3[2], tf3[3], tf3[4])
00062 tf=geom.Mat4(rot)
00063 tf.PasteTranslation(geom.Vec3(tf1[1], tf2[1], tf3[1]))
00064 return TMAlignResult(rmsd, aln_length, tm_score, tf)
00065
00066 def _RunTmAlign(tmalign, tmp_dir):
00067 model1_filename=os.path.join(tmp_dir, 'model01.pdb')
00068 model2_filename=os.path.join(tmp_dir, 'model02.pdb')
00069 if platform.system() == "Windows":
00070 tmalign_path=settings.Locate('tmalign.exe', explicit_file_name=tmalign)
00071 command="\"%s\" %s %s" %(os.path.normpath(tmalign_path), model1_filename, model2_filename)
00072 else:
00073 tmalign_path=settings.Locate('tmalign', explicit_file_name=tmalign)
00074 command="%s '%s' '%s'" %(tmalign_path, model1_filename, model2_filename)
00075 ps=subprocess.Popen(command, shell=True, stdout=subprocess.PIPE)
00076 ps.wait()
00077 lines=ps.stdout.readlines()
00078 if (len(lines))<22:
00079 raise RuntimeError("tmalign superposition failed")
00080 return _ParseTmAlign(lines)
00081
00082 class TMScoreResult:
00083 def __init__(self, rmsd_common, tm_score, max_sub,
00084 gdt_ts, gdt_ha, rmsd_below_five, transform):
00085 self.rmsd_common=rmsd_common
00086 self.tm_score=tm_score
00087 self.max_sub=max_sub
00088 self.gdt_ts=gdt_ts
00089 self.gdt_ha=gdt_ha
00090 self.rmsd_below_five=rmsd_below_five
00091 self.transform=transform
00092
00093 def _ParseTmScore(lines):
00094 tf1=[float(i.strip()) for i in lines[23].split()]
00095 tf2=[float(i.strip()) for i in lines[24].split()]
00096 tf3=[float(i.strip()) for i in lines[25].split()]
00097 rot=geom.Mat3(tf1[2], tf1[3], tf1[4], tf2[2], tf2[3],
00098 tf2[4], tf3[2], tf3[3], tf3[4])
00099 tf=geom.Mat4(rot)
00100 tf.PasteTranslation(geom.Vec3(tf1[1], tf2[1], tf3[1]))
00101 result=TMScoreResult(float(lines[14].split()[-1].strip()),
00102 float(lines[16].split()[2].strip()),
00103 float(lines[17].split()[1].strip()),
00104 float(lines[18].split()[1].strip()),
00105 float(lines[19].split()[1].strip()),
00106 float(lines[27].split()[-1].strip()),
00107 tf)
00108 return result
00109
00110
00111 def _RunTmScore(tmscore, tmp_dir):
00112 model1_filename=os.path.join(tmp_dir, 'model01.pdb')
00113 model2_filename=os.path.join(tmp_dir, 'model02.pdb')
00114 if platform.system() == "Windows":
00115 tmscore_path=settings.Locate('tmscore.exe', explicit_file_name=tmscore)
00116 command="\"%s\" %s %s" %(os.path.normpath(tmscore_path), model1_filename, model2_filename)
00117 else:
00118 tmscore_path=settings.Locate('tmscore', explicit_file_name=tmscore)
00119 command="%s '%s' '%s'" %(tmscore_path, model1_filename, model2_filename)
00120 ps=subprocess.Popen(command, shell=True, stdout=subprocess.PIPE)
00121 ps.wait()
00122 lines=ps.stdout.readlines()
00123 if (len(lines))<22:
00124 raise RuntimeError("tmscore superposition failed")
00125 return _ParseTmScore(lines)
00126
00127 def TMAlign(model1, model2, tmalign=None):
00128 """
00129 Run tmalign on two protein structures
00130 """
00131 tmp_dir_name=_SetupFiles((model1, model2))
00132 result=_RunTmAlign(tmalign, tmp_dir_name)
00133 model1.handle.RequestXCSEditor().ApplyTransform(result.transform)
00134 _CleanupFiles(tmp_dir_name)
00135 return result
00136 def TMScore(model1, model2, tmscore=None):
00137 """
00138 Run tmscore on two protein structures
00139 """
00140 tmp_dir_name=_SetupFiles((model1, model2))
00141 result=_RunTmScore(tmscore, tmp_dir_name)
00142 model1.handle.RequestXCSEditor().ApplyTransform(result.transform)
00143 _CleanupFiles(tmp_dir_name)
00144 return result