mirror of
https://github.com/schrodinger/pymol-open-source.git
synced 2026-06-03 19:54:24 +08:00
Nucleic Acid builder and general improvements (#354)
Fixes #102 #351 Migrating support from Incentive for the following: - Nucleic Acid builder in the Builder panel (Jarrett Johnson @JarrettSJohnson) - `fnab` command to build nucleic acid chains by sequence (Jarrett Johnson @JarrettSJohnson) - Improvements to nucleic acid building from native structures (Thomas Stewart @TstewDev) - highlight attachment points for improved usability (Thomas Holder @speleo3) Co-authored-by: Jarrett Johnson <jarrett.johnson@schrodinger.com> Co-authored-by: Thomas Stewart <thomas.stewart@schrodinger.com> Co-authored-by: Thomas Holder <thomas@thomas-holder.de>
This commit is contained in:
BIN
data/chempy/fragments/atpA.pkl
Normal file
BIN
data/chempy/fragments/atpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/atpB.pkl
Normal file
BIN
data/chempy/fragments/atpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/atp_ttpA.pkl
Normal file
BIN
data/chempy/fragments/atp_ttpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/atp_ttpB.pkl
Normal file
BIN
data/chempy/fragments/atp_ttpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ctpA.pkl
Normal file
BIN
data/chempy/fragments/ctpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ctpB.pkl
Normal file
BIN
data/chempy/fragments/ctpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ctp_gtpA.pkl
Normal file
BIN
data/chempy/fragments/ctp_gtpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ctp_gtpB.pkl
Normal file
BIN
data/chempy/fragments/ctp_gtpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/gtpA.pkl
Normal file
BIN
data/chempy/fragments/gtpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/gtpB.pkl
Normal file
BIN
data/chempy/fragments/gtpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/gtp_ctpA.pkl
Normal file
BIN
data/chempy/fragments/gtp_ctpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/gtp_ctpB.pkl
Normal file
BIN
data/chempy/fragments/gtp_ctpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ttpA.pkl
Normal file
BIN
data/chempy/fragments/ttpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ttpB.pkl
Normal file
BIN
data/chempy/fragments/ttpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ttp_atpA.pkl
Normal file
BIN
data/chempy/fragments/ttp_atpA.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/ttp_atpB.pkl
Normal file
BIN
data/chempy/fragments/ttp_atpB.pkl
Normal file
Binary file not shown.
BIN
data/chempy/fragments/utpA.pkl
Normal file
BIN
data/chempy/fragments/utpA.pkl
Normal file
Binary file not shown.
@@ -13,7 +13,6 @@ from pymol.Qt.utils import PopupOnException
|
||||
Qt = QtCore.Qt
|
||||
|
||||
from pymol.wizard import Wizard
|
||||
from pymol.parsing import QuietException
|
||||
|
||||
from pymol import editor
|
||||
from pymol import computing
|
||||
@@ -306,7 +305,7 @@ class AttachWizard(RepeatableActionWizard):
|
||||
RepeatableActionWizard.__init__(self,_self)
|
||||
self.mode = 0
|
||||
|
||||
def do_pick(self, bondFlag):
|
||||
def do_pick(self, bondFlag,*,_self=pymol.cmd):
|
||||
if self.mode == 0:
|
||||
self.cmd.select(active_sele, "bymol pk1")
|
||||
editor.attach_fragment("pk1", self.fragment, self.position, self.geometry, _self=self.cmd)
|
||||
@@ -328,14 +327,9 @@ class AttachWizard(RepeatableActionWizard):
|
||||
self.setActionHash( (fragment, position, geometry, text) )
|
||||
self.activateRepeatOrDismiss()
|
||||
|
||||
def create_new(self):
|
||||
names = self.cmd.get_names("objects")
|
||||
num = 1
|
||||
while 1:
|
||||
name = "obj%02d"%num
|
||||
if name not in names:
|
||||
break
|
||||
num = num + 1
|
||||
def create_new(self,*,_self=pymol.cmd):
|
||||
self.cmd.unpick()
|
||||
name = self.cmd.get_unused_name("obj")
|
||||
self.cmd.fragment(self.fragment, name)
|
||||
if not self.getRepeating():
|
||||
self.actionWizardDone()
|
||||
@@ -356,7 +350,7 @@ class AttachWizard(RepeatableActionWizard):
|
||||
def get_panel(self):
|
||||
if self.getRepeating():
|
||||
return [
|
||||
[ 1, 'Attaching Multiple Fragmnts',''],
|
||||
[ 1, 'Attaching Multiple Fragments',''],
|
||||
[ 2, 'Create As New Object','cmd.get_wizard().create_new()'],
|
||||
[ 2, 'Combine w/ Existing Object','cmd.get_wizard().combine()'],
|
||||
[ 2, 'Done','cmd.set_wizard()'],
|
||||
@@ -371,34 +365,56 @@ class AttachWizard(RepeatableActionWizard):
|
||||
]
|
||||
|
||||
|
||||
class AminoAcidWizard(RepeatableActionWizard):
|
||||
class BioPolymerWizard(RepeatableActionWizard):
|
||||
|
||||
def __init__(self, _self=pymol.cmd, ss=-1):
|
||||
HIGHLIGHT_SELE = ""
|
||||
|
||||
def __init__(self, _self=pymol.cmd):
|
||||
RepeatableActionWizard.__init__(self,_self)
|
||||
self.mode = 0
|
||||
self.setSecondaryStructure(ss)
|
||||
self._highlighting_enabled = False
|
||||
|
||||
def setSecondaryStructure(self, ss):
|
||||
self._secondary_structure = ss
|
||||
def __enter__(self):
|
||||
'''Context manager for temporarily disabling attachment point highlights
|
||||
'''
|
||||
self.highlight_attachment_points(False)
|
||||
|
||||
def __exit__(self, *_):
|
||||
self.highlight_attachment_points()
|
||||
|
||||
def cleanup(self):
|
||||
self.highlight_attachment_points(False)
|
||||
RepeatableActionWizard.cleanup(self)
|
||||
|
||||
def highlight_attachment_points(self, show=True):
|
||||
'''Show spheres for potential attachment points
|
||||
|
||||
:type show: bool
|
||||
:param show: Switch to show or hide the highlights
|
||||
'''
|
||||
if self._highlighting_enabled:
|
||||
fn = self.cmd.show if show else self.cmd.hide
|
||||
fn('spheres', self.HIGHLIGHT_SELE)
|
||||
|
||||
def attach_monomer(self, objectname=""):
|
||||
editor.attach_amino_acid("?pk1", self.aminoAcid, object=objectname,
|
||||
ss=self._secondary_structure,
|
||||
_self=self.cmd)
|
||||
raise NotImplementedError
|
||||
|
||||
def combine_monomer(self):
|
||||
raise NotImplementedError
|
||||
|
||||
def do_pick(self, bondFlag):
|
||||
# since this function can change any position of atoms in a related
|
||||
# molecule, bymol is used
|
||||
if self.mode == 0:
|
||||
self.cmd.select(active_sele, "bymol pk1")
|
||||
self.cmd.select(active_sele, "bymol ?pk1")
|
||||
try:
|
||||
with undocontext(self.cmd, "bymol ?pk1"):
|
||||
self.attach_monomer(self.aminoAcid)
|
||||
except QuietException:
|
||||
fin = -1
|
||||
self.attach_monomer()
|
||||
except pymol.CmdException as exc:
|
||||
print(exc)
|
||||
elif self.mode == 1:
|
||||
self.cmd.select(active_sele, "bymol pk1")
|
||||
editor.combine_fragment("pk1", self.aminoAcid, 0, 1, _self=self.cmd)
|
||||
self.cmd.select(active_sele, "bymol ?pk1")
|
||||
editor.combine_monomer()
|
||||
self.mode = 0
|
||||
self.cmd.refresh_wizard()
|
||||
|
||||
@@ -406,23 +422,25 @@ class AminoAcidWizard(RepeatableActionWizard):
|
||||
if not self.getRepeating():
|
||||
self.actionWizardDone()
|
||||
|
||||
def toggle(self, amino_acid):
|
||||
self.aminoAcid = amino_acid
|
||||
self.setActionHash( (amino_acid,) )
|
||||
self.activateRepeatOrDismiss()
|
||||
def toggle(self, monomer):
|
||||
self._monomer = monomer
|
||||
self.setActionHash( (monomer,) )
|
||||
if self.activateRepeatOrDismiss():
|
||||
# enable temporary sphere highlighting of attachment points, but only
|
||||
# if currently no spheres are displayed
|
||||
self._highlighting_enabled = self.HIGHLIGHT_SELE and self.cmd.count_atoms(
|
||||
'(rep spheres) & ({})'.format(self.HIGHLIGHT_SELE)) == 0
|
||||
|
||||
def create_new(self):
|
||||
names = self.cmd.get_names("objects")
|
||||
num = 1
|
||||
while 1:
|
||||
name = "obj%02d"%num
|
||||
if name not in names:
|
||||
break
|
||||
num = num + 1
|
||||
self.attach_monomer(self.aminoAcid)
|
||||
self.highlight_attachment_points()
|
||||
|
||||
def create_new(self, *, _self=pymol.cmd):
|
||||
self.cmd.unpick()
|
||||
name = self.cmd.get_unused_name("obj")
|
||||
self.attach_monomer(name)
|
||||
if not self.getRepeating():
|
||||
self.actionWizardDone()
|
||||
else:
|
||||
self.cmd.unpick()
|
||||
|
||||
def combine(self):
|
||||
self.mode = 1
|
||||
@@ -431,29 +449,67 @@ class AminoAcidWizard(RepeatableActionWizard):
|
||||
def get_prompt(self):
|
||||
if self.mode == 0:
|
||||
if self.getRepeating():
|
||||
return ["Pick locations to attach %s..."%self.aminoAcid]
|
||||
return ["Pick locations to attach %s..." % self._monomer]
|
||||
else:
|
||||
return ["Pick location to attach %s..."%self.aminoAcid]
|
||||
return ["Pick location to attach %s..." % self._monomer]
|
||||
else:
|
||||
return ["Pick object to combine %s into..."%self.aminoAcid]
|
||||
return ["Pick object to combine %s into..." % self._monomer]
|
||||
|
||||
def get_panel(self):
|
||||
if self.getRepeating():
|
||||
return [
|
||||
[ 1, 'Attaching Multiple Residues',''],
|
||||
[ 2, 'Create As New Object','cmd.get_wizard().create_new()'],
|
||||
[ 2, 'Combine w/ Existing Object','cmd.get_wizard().combine()'],
|
||||
[ 2, 'Done','cmd.set_wizard()'],
|
||||
]
|
||||
else:
|
||||
return [
|
||||
[ 1, 'Attaching Amino Acid',''],
|
||||
[ 2, 'Create As New Object','cmd.get_wizard().create_new()'],
|
||||
[ 2, 'Combine w/ Existing Object','cmd.get_wizard().combine()'],
|
||||
[ 2, 'Attach Multiple...','cmd.get_wizard().repeat()'],
|
||||
[ 2, 'Done','cmd.set_wizard()'],
|
||||
]
|
||||
|
||||
class AminoAcidWizard(BioPolymerWizard):
|
||||
|
||||
HIGHLIGHT_SELE = "(name N &! neighbor name C) | (name C &! neighbor name N)"
|
||||
|
||||
def __init__(self, _self=pymol.cmd, ss=-1):
|
||||
BioPolymerWizard.__init__(self,_self)
|
||||
self._monomerType = "Amino Acid"
|
||||
self.setSecondaryStructure(ss)
|
||||
|
||||
def setSecondaryStructure(self, ss):
|
||||
self._secondary_structure = ss
|
||||
|
||||
def attach_monomer(self, objectname=""):
|
||||
with self:
|
||||
editor.attach_amino_acid("?pk1", self._monomer, object=objectname,
|
||||
ss=self._secondary_structure,
|
||||
_self=self.cmd)
|
||||
|
||||
def combine_monomer(self):
|
||||
editor.combine_fragment("pk1", self._monomer, 0, 1, _self=self.cmd)
|
||||
|
||||
class NucleicAcidWizard(BioPolymerWizard):
|
||||
|
||||
HIGHLIGHT_SELE = "(name O3' &! neighbor name P) | (name P &! neighbor name O3') | (name O5' &! neighbor name P) "
|
||||
|
||||
def _init(self, form, dbl_helix, nuc_type):
|
||||
self._monomerType = "Nucleic Acid"
|
||||
self._form = form
|
||||
self._dbl_helix = dbl_helix
|
||||
self._nuc_type = nuc_type
|
||||
return self
|
||||
|
||||
def attach_monomer(self, objectname=""):
|
||||
with self:
|
||||
editor.attach_nuc_acid("?pk1", self._monomer, object=objectname,
|
||||
nuc_type=self._nuc_type, form=self._form,
|
||||
dbl_helix=self._dbl_helix, _self=self.cmd)
|
||||
|
||||
def combine_monomer(self):
|
||||
editor.combine_nucleotide("pk1", self._monomer + self._form, 0, 1, _self=self.cmd)
|
||||
|
||||
class ValenceWizard(RepeatableActionWizard):
|
||||
|
||||
@@ -950,6 +1006,14 @@ def collectPicked(self_cmd):
|
||||
return result
|
||||
|
||||
|
||||
#############################################################
|
||||
### Nucleic Acid helper class
|
||||
|
||||
class NucleicAcidProperties:
|
||||
def __init__(self, form="B", double_helix=True):
|
||||
self.dna_form = form
|
||||
self.dna_dbl_helix = double_helix
|
||||
|
||||
#############################################################
|
||||
### Actual GUI
|
||||
|
||||
@@ -995,8 +1059,15 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
self.protein_tab = QtWidgets.QWidget()
|
||||
self.protein_tab.setLayout(self.protein_layout)
|
||||
|
||||
self.nucleic_acid_layout = QtWidgets.QVBoxLayout()
|
||||
self.nucleic_acid_layout.setContentsMargins(5, 5, 5, 5)
|
||||
self.nucleic_acid_layout.setSpacing(5)
|
||||
self.nucleic_acid_tab = QtWidgets.QWidget()
|
||||
self.nucleic_acid_tab.setLayout(self.nucleic_acid_layout)
|
||||
|
||||
self.tabs.addTab(self.fragments_tab, "Chemical")
|
||||
self.tabs.addTab(self.protein_tab, "Protein")
|
||||
self.tabs.addTab(self.nucleic_acid_tab, "Nucleic Acid")
|
||||
|
||||
self.getIcons()
|
||||
|
||||
@@ -1082,6 +1153,76 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
self.protein_layout.addWidget(self.ss_cbox, 2, lab_cols, 1, 4)
|
||||
self.ss_cbox.currentIndexChanged[int].connect(self.ssIndexChanged)
|
||||
|
||||
self.nucleic_acid_dna_layout = QtWidgets.QGridLayout()
|
||||
self.nucleic_acid_dna_layout.setContentsMargins(5, 5, 5, 5)
|
||||
self.nucleic_acid_dna_layout.setSpacing(5)
|
||||
self.nucleic_acid_rna_layout = QtWidgets.QGridLayout()
|
||||
self.nucleic_acid_rna_layout.setContentsMargins(5, 5, 5, 5)
|
||||
self.nucleic_acid_rna_layout.setSpacing(5)
|
||||
|
||||
|
||||
self.nucleic_acid_tab = QtWidgets.QTabWidget()
|
||||
self.nucleic_acid_layout.addWidget(self.nucleic_acid_tab)
|
||||
self.dna_tab = QtWidgets.QWidget()
|
||||
self.dna_tab.setLayout(self.nucleic_acid_dna_layout)
|
||||
self.rna_tab = QtWidgets.QWidget()
|
||||
self.rna_tab.setLayout(self.nucleic_acid_rna_layout)
|
||||
|
||||
self.nucleic_acid_tab.addTab(self.dna_tab, "DNA")
|
||||
self.nucleic_acid_tab.addTab(self.rna_tab, "RNA")
|
||||
|
||||
self._nuc_acid_prop = NucleicAcidProperties()
|
||||
|
||||
dna_buttons = (
|
||||
("A", "Deoxyadenosine", None, lambda: self.attach_nuc_acid("atp", "DNA")),
|
||||
("C", "Deoxycytidine", None, lambda: self.attach_nuc_acid("ctp", "DNA")),
|
||||
("T", "Deoxythymidine", None, lambda: self.attach_nuc_acid("ttp", "DNA")),
|
||||
("G", "Deoxyguanosine", None, lambda: self.attach_nuc_acid("gtp", "DNA")),
|
||||
("@Form:", None, None, None),
|
||||
("#A", None, False, lambda: setattr(self._nuc_acid_prop, 'dna_form', "A")),
|
||||
("#B", None, True, lambda: setattr(self._nuc_acid_prop, 'dna_form', "B")),
|
||||
("@Helix:", None, None, None),
|
||||
("#Single", None, False, lambda: setattr(self._nuc_acid_prop, 'dna_dbl_helix', False)),
|
||||
("#Double", None, True, lambda: setattr(self._nuc_acid_prop, 'dna_dbl_helix', True))
|
||||
)
|
||||
for col_num, btn_pkg in enumerate(dna_buttons):
|
||||
btn_label, btn_tooltip, default_activated, btn_command = btn_pkg
|
||||
if btn_label[0] == '@':
|
||||
btn = QtWidgets.QLabel(btn_label[1:])
|
||||
radio_group = QtWidgets.QButtonGroup(self)
|
||||
elif btn_label[0] == '#':
|
||||
btn = QtWidgets.QRadioButton(btn_label[1:])
|
||||
btn.toggled.connect(btn_command)
|
||||
btn.setChecked(default_activated)
|
||||
radio_group.addButton(btn)
|
||||
else:
|
||||
btn = makeFragmentButton()
|
||||
btn.setText(btn_label)
|
||||
btn.setToolTip(btn_tooltip)
|
||||
btn.clicked.connect(btn_command)
|
||||
self.nucleic_acid_dna_layout.addWidget(btn, 0, col_num)
|
||||
|
||||
|
||||
rna_buttons = (
|
||||
("A", "Adenosine", "atp", lambda: self.attach_nuc_acid("atp", "RNA")),
|
||||
("C", "Cytosine", "ctp", lambda: self.attach_nuc_acid("ctp", "RNA")),
|
||||
("U", "Uracil", "utp", lambda: self.attach_nuc_acid("utp", "RNA")),
|
||||
("G", "Guanine", "gtp", lambda: self.attach_nuc_acid("gtp", "RNA")))
|
||||
|
||||
for col_num, btn_pkg in enumerate(rna_buttons):
|
||||
btn_label, btn_tooltip, btn_filename, btn_command = btn_pkg
|
||||
btn = makeFragmentButton()
|
||||
btn.setText(btn_label)
|
||||
btn.setToolTip(btn_tooltip)
|
||||
btn.clicked.connect(btn_command)
|
||||
self.nucleic_acid_rna_layout.addWidget(btn, 0, col_num)
|
||||
|
||||
btn = QtWidgets.QLabel('Hint: Also check out '
|
||||
'<a href="http://x3dna.org/articles/3dna-fiber-models">fiber</a> and its '
|
||||
'<a href="http://x3dna.org/articles/pymol-wrapper-to-3dna-fiber-models">PyMOL wrapper</a>')
|
||||
btn.setOpenExternalLinks(True)
|
||||
self.nucleic_acid_rna_layout.addWidget(btn, 0, col_num + 1)
|
||||
|
||||
buttons = [
|
||||
[
|
||||
( "@Atoms:", None, None),
|
||||
@@ -1094,6 +1235,8 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
( " +1 ", "Positive Charge", lambda: self.setCharge(1,"+1")),
|
||||
( " 0 ", "Neutral Charge", lambda: self.setCharge(0,"neutral")),
|
||||
( " -1 ", "Negative Charge", lambda: self.setCharge(-1,"-1")),
|
||||
( "@ Residue:", None, None),
|
||||
("Remove", "Remove residue", lambda: self.removeResn()),
|
||||
],
|
||||
[
|
||||
( "@Bonds:", None, None),
|
||||
@@ -1139,9 +1282,12 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
self.cmd.set(n, checked, quiet=0)
|
||||
else:
|
||||
btn.setChecked(not value)
|
||||
@btn.toggled.connect
|
||||
def _(checked, n=setting):
|
||||
self.cmd.set(n, not checked, quiet=0)
|
||||
if setting == 'suspend_undo':
|
||||
self.setUndoEnabled(not value)
|
||||
_ = self.setUndoEnabled
|
||||
btn.toggled.connect(_)
|
||||
else:
|
||||
btn = makeFragmentButton()
|
||||
btn.setText(btn_label)
|
||||
@@ -1151,6 +1297,29 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
btn_row_layout.addWidget(btn)
|
||||
btn_row_layout.addStretch()
|
||||
|
||||
def setUndoEnabled(self, checked):
|
||||
self.cmd.set('suspend_undo', not checked, quiet=0)
|
||||
if not checked:
|
||||
return
|
||||
|
||||
on_per_object = set(oname for oname in self.cmd.get_object_list()
|
||||
if self.cmd.get_setting_int('suspend_undo', oname))
|
||||
|
||||
n = len(on_per_object)
|
||||
if n > 20:
|
||||
on_per_object = sorted(on_per_object)[:15] + [
|
||||
"[{} more]".format(n - 15)]
|
||||
|
||||
if on_per_object:
|
||||
QMB = QtWidgets.QMessageBox
|
||||
check = QMB.question(None, 'Enable for objects?',
|
||||
'Building "Undo" is disabled for the following objects:\n\n' +
|
||||
'\n'.join(on_per_object) + '\n\n'
|
||||
'Enable "Undo" for these objects?', QMB.Yes | QMB.No)
|
||||
if check == QMB.Yes:
|
||||
for oname in on_per_object:
|
||||
self.cmd.unset('suspend_undo', oname)
|
||||
|
||||
def getIcons(self):
|
||||
self.icons = {}
|
||||
# use old Tk icons
|
||||
@@ -1209,6 +1378,37 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
if isinstance(w, AminoAcidWizard):
|
||||
w.setSecondaryStructure(index + 1)
|
||||
|
||||
def attach_nuc_acid(self, nuc_acid, nuc_type):
|
||||
self._nuc_type = nuc_type
|
||||
picked = collectPicked(self.cmd)
|
||||
if len(picked)==1:
|
||||
try:
|
||||
with undocontext(self.cmd, "byobject %s" % picked[0]):
|
||||
editor.attach_nuc_acid(picked[0], nuc_acid,
|
||||
nuc_type=self._nuc_type,
|
||||
object="",
|
||||
form=self._nuc_acid_prop.dna_form,
|
||||
dbl_helix=self._nuc_acid_prop.dna_dbl_helix,
|
||||
_self=self.cmd)
|
||||
except pymol.CmdException as exc:
|
||||
print(exc)
|
||||
except ValueError as exc:
|
||||
print(exc)
|
||||
self.doZoom()
|
||||
else:
|
||||
self.cmd.unpick()
|
||||
NucleicAcidWizard(_self=self.cmd)._init(form=self._nuc_acid_prop.dna_form,
|
||||
dbl_helix=self._nuc_acid_prop.dna_dbl_helix,
|
||||
nuc_type=self._nuc_type).toggle(nuc_acid)
|
||||
|
||||
def removeResn(self):
|
||||
picked = collectPicked(self.cmd)
|
||||
if picked == ["pk1"]:
|
||||
self.cmd.select(newest_sele,"byres(pk1)")
|
||||
self.cmd.remove(newest_sele)
|
||||
else:
|
||||
print("Select a single atom on the residue and press remove again")
|
||||
|
||||
def doAutoPick(self, old_atoms=None):
|
||||
self.cmd.unpick()
|
||||
if self.cmd.select(newest_sele,"(byobj "+active_sele+") and not "+active_sele)==0:
|
||||
@@ -1228,7 +1428,7 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
|
||||
def doZoom(self, *ignore):
|
||||
if "pk1" in self.cmd.get_names("selections"):
|
||||
self.cmd.zoom("((neighbor pk1) extend 4)", 4.0, animate=-1)
|
||||
self.cmd.center("%pk1 extend 9", animate=-1)
|
||||
|
||||
def setCharge(self, charge, text):
|
||||
picked = collectPicked(self.cmd)
|
||||
@@ -1304,12 +1504,6 @@ class _BuilderPanel(QtWidgets.QWidget):
|
||||
self.cmd.unpick()
|
||||
InvertWizard(self.cmd).toggle()
|
||||
|
||||
def center(self):
|
||||
if "pk1" in self.cmd.get_names("selections"):
|
||||
self.cmd.zoom("pk1", 5.0, animate=-1)
|
||||
else:
|
||||
self.cmd.zoom("all", 3.0, animate=-1)
|
||||
|
||||
def removeAtom(self):
|
||||
picked = collectPicked(self.cmd)
|
||||
if len(picked):
|
||||
|
||||
@@ -345,7 +345,7 @@ class AttachWizard(RepeatableActionWizard):
|
||||
def get_panel(self):
|
||||
if self.getRepeating():
|
||||
return [
|
||||
[ 1, 'Attaching Multiple Fragmnts',''],
|
||||
[ 1, 'Attaching Multiple Fragments',''],
|
||||
[ 2, 'Create As New Object','cmd.get_wizard().create_new()'],
|
||||
[ 2, 'Combine w/ Existing Object','cmd.get_wizard().combine()'],
|
||||
[ 2, 'Done','cmd.set_wizard()'],
|
||||
|
||||
@@ -259,7 +259,8 @@ from .editing import \
|
||||
vdw_fit
|
||||
|
||||
from .editor import \
|
||||
fab
|
||||
fab, \
|
||||
fnab
|
||||
|
||||
from .computing import \
|
||||
clean
|
||||
|
||||
@@ -1,3 +1,5 @@
|
||||
import math
|
||||
|
||||
import re
|
||||
import pymol
|
||||
cmd = __import__("sys").modules["pymol.cmd"]
|
||||
@@ -384,6 +386,659 @@ def _fab(input,name,mode,resi,chain,segi,state,dir,hydro,ss,quiet,_self=cmd):
|
||||
|
||||
return r
|
||||
|
||||
_threeNA_to_OneNA = { "atp" : "A", "ctp" : "C", "gtp" : "G", "ttp" : "T", "utp" : "U"}
|
||||
_oneNA_to_threeNA = { "A" : "atp", "C" : "ctp", "G" : "gtp", "T" : "ttp", "U" : "utp"}
|
||||
|
||||
_base_pair = { "DNA" : {"atp" : "ttp", "ctp" : "gtp", "gtp" : "ctp", "ttp" : "atp" },
|
||||
"RNA" : {"atp" : "utp", "ctp" : "gtp", "gtp" : "ctp", "utp" : "atp" }}
|
||||
|
||||
_oneNA_base_pair = {"DA" : "DT", "DC" : "DG", "DG" : "DC", "DT" : "DA", "A" : "U",
|
||||
"C" : "G", "G" : "C", "U" : "A"}
|
||||
|
||||
def iterate_to_list(selection: str, expression: str, *, _self=cmd):
|
||||
outlist = []
|
||||
_self.iterate(selection,f"outlist.append(({expression}))", space={"outlist":outlist})
|
||||
return outlist
|
||||
|
||||
def rename_three_to_one(nuc_acid, sele, nuc_type, *, _self=cmd):
|
||||
"""
|
||||
Renames nucleobase from 3-letter to 1-letter representation
|
||||
|
||||
:param nuc_acid: (str) 3-letter nucleic acid representation
|
||||
:param sele: (str) selection of nucleic acid to rename
|
||||
:param nuc_type: (str) "DNA" or "RNA"
|
||||
"""
|
||||
new_name = _threeNA_to_OneNA[nuc_acid]
|
||||
if nuc_type == "DNA":
|
||||
new_name = "D" + new_name
|
||||
_self.alter(sele, f"resn='{new_name}'")
|
||||
|
||||
def fit_sugars(mobile, target, *, _self=cmd):
|
||||
"""
|
||||
Fits appending base pairs to form appropriate hydrogen bond
|
||||
|
||||
:param mobile: (str) selection for the sense (main) strand
|
||||
:param target: (str) selection for the antisense (opposing) strand
|
||||
"""
|
||||
try:
|
||||
_self.pair_fit(f"{mobile} & name C1'",
|
||||
f"{target} & name C1'",
|
||||
f"{mobile} & name C2'",
|
||||
f"{target} & name C2'",
|
||||
f"{mobile} & name C3'",
|
||||
f"{target} & name C3'",
|
||||
f"{mobile} & name C4'",
|
||||
f"{target} & name C4'",
|
||||
f"{mobile} & name O4'",
|
||||
f"{target} & name O4'", quiet=1)
|
||||
except:
|
||||
_self.delete(tmp_wild)
|
||||
raise pymol.CmdException("Something went wrong when fitting the new residue.")
|
||||
|
||||
def fit_DS_fragment(mobile_A, target_A, mobile_B, target_B, *, _self=cmd):
|
||||
"""
|
||||
Fits dummy fragment to the detected structure using atoms
|
||||
on both stands for a more accurate alignment.
|
||||
|
||||
:param mobile_A: (str) selection for the base being created and attached
|
||||
:param target_A: (str) selection for the base selected to build on
|
||||
:param mobile_B: (str) selection for the opposing base being created
|
||||
:param target_B: (str) selection for the detected opposing base
|
||||
"""
|
||||
try:
|
||||
_self.pair_fit(f"{mobile_A} & name C1'",
|
||||
f"{target_A} & name C1'",
|
||||
f"{mobile_A} & name C2'",
|
||||
f"{target_A} & name C2'",
|
||||
f"{mobile_A} & name C5'",
|
||||
f"{target_A} & name C5'",
|
||||
f"{mobile_A} & name O4'",
|
||||
f"{target_A} & name O4'",
|
||||
f"{mobile_A} & name O3'",
|
||||
f"{target_A} & name O3'",
|
||||
f"{mobile_A} & name P",
|
||||
f"{target_A} & name P",
|
||||
f"{mobile_B} & name O3'",
|
||||
f"{target_B} & name O3'",
|
||||
f"{mobile_B} & name C1'",
|
||||
f"{target_B} & name C1'",
|
||||
f"{mobile_B} & name C2'",
|
||||
f"{target_B} & name C2'",
|
||||
f"{mobile_B} & name C5'",
|
||||
f"{target_B} & name C5'",
|
||||
f"{mobile_B} & name P",
|
||||
f"{target_B} & name P",
|
||||
f"{mobile_B} & name O4'",
|
||||
f"{target_B} & name O4'", quiet=1)
|
||||
except:
|
||||
_self.delete(tmp_wild)
|
||||
raise pymol.CmdException("Something went wrong when fitting the new residue.")
|
||||
|
||||
def add2pO(domain, nuc_acid, resv, *, _self=cmd):
|
||||
if nuc_acid == "utp": #utp comes with O2'
|
||||
return
|
||||
c_2p = _prefix + "_c2p"
|
||||
_self.select(c_2p, "%s & resi \\%i & name %s" % (domain, resv, "C2'"))
|
||||
_self.unpick()
|
||||
_self.edit(c_2p)
|
||||
_self.attach("O", 4, 4)
|
||||
_self.unpick()
|
||||
_self.alter("(byres %s) & resi \\%i & name O01" % (domain, resv), "name=\"O2'\"")
|
||||
|
||||
def move_atom_in_res(atom, dummy_res, new_res, twist, rise, *, _self=cmd):
|
||||
prev_coords = _self.get_coords(f"{dummy_res} & name {atom}", state=1)
|
||||
curr_coords = _self.get_coords(f"{new_res} & name {atom}", state=1)
|
||||
|
||||
if curr_coords is None or prev_coords is None:
|
||||
return
|
||||
curr_coord = curr_coords[0]
|
||||
prev_coord = prev_coords[0]
|
||||
|
||||
r = math.sqrt(prev_coord[0] ** 2 +
|
||||
prev_coord[1] ** 2)
|
||||
old_phi = math.degrees(math.atan2(prev_coord[1],prev_coord[0]))
|
||||
phi = old_phi - twist
|
||||
phi = math.radians(phi)
|
||||
new_pos = [r * math.cos(phi),
|
||||
r * math.sin(phi),
|
||||
prev_coord[2] - rise]
|
||||
|
||||
trans = list(new_pos - curr_coord)
|
||||
_self.translate(trans, f"{new_res} & name {atom}", camera=0)
|
||||
|
||||
def move_new_res(frag_string, full_frag, old, old_oppo, double_stranded_bool=False, form="B", chain="A", antisense=False, *, _self=cmd):
|
||||
"""
|
||||
Attaches new residue (or pair) onto current nucleotide chain
|
||||
|
||||
:param frag_string: (str) Name of appending nucleic acid or pair
|
||||
:param full_frag: (str) Selection of the created fragment
|
||||
:param old: (str) Selection of previous residue
|
||||
:param old_oppo: (str) Selection of previous opposing residue
|
||||
:param double_stranded_bool: (bool) Flag represing if double helix was detected
|
||||
:param antisense (bool) Flag for antisense
|
||||
:param form: (str) DNA form ('A'/'B')
|
||||
"""
|
||||
if form == 'B':
|
||||
twist = -36.0
|
||||
rise = -3.375
|
||||
elif form == 'A':
|
||||
twist = -32.7
|
||||
rise = -2.548
|
||||
else:
|
||||
raise ValueError("Form not recognized")
|
||||
|
||||
rise = rise if antisense else -rise
|
||||
twist = twist if antisense else -twist
|
||||
|
||||
dummy_fragment = _prefix + "_dummyfrag"
|
||||
dummy_res_A = _prefix + "_dummyresA"
|
||||
dummy_res_B = _prefix + "_dummyresB"
|
||||
new_fragment = _prefix + "_newfrag"
|
||||
new_res_A = _prefix + "_newresA"
|
||||
new_res_B = _prefix + "_newresB"
|
||||
|
||||
_self.select(new_fragment, f"{full_frag}")
|
||||
_self.select(new_res_A,f"{full_frag} and chain A")
|
||||
_self.select(new_res_B,f"{full_frag} and chain B")
|
||||
|
||||
_self.fragment(frag_string, dummy_fragment, origin=0)
|
||||
_self.select(dummy_res_A, f"{dummy_fragment} and chain A")
|
||||
_self.select(dummy_res_B, f"{dummy_fragment} and chain B")
|
||||
|
||||
if old_oppo == "none":
|
||||
# This is the case where a single residue is being added
|
||||
_self.select(dummy_res_A, f"{dummy_fragment}")
|
||||
_self.select(new_res_A, f"{full_frag}")
|
||||
|
||||
atoms_A = iterate_to_list(dummy_res_A, "name")
|
||||
atoms_B = iterate_to_list(dummy_res_B, "name")
|
||||
|
||||
#A new base is created by copying the coordinates of the previous
|
||||
#base and doing a cylindrical rotation (phi degrees) and a translation
|
||||
#down the z-axis by the rise amount
|
||||
for atom in atoms_A:
|
||||
move_atom_in_res(atom, dummy_res_A, new_res_A, twist, rise)
|
||||
for atom in atoms_B:
|
||||
move_atom_in_res(atom, dummy_res_B, new_res_B, twist, rise)
|
||||
if double_stranded_bool == True:
|
||||
fit_DS_fragment(dummy_res_A, old, dummy_res_B, old_oppo)
|
||||
orient_flag = check_dummy_oriention(old,dummy_res_A)
|
||||
if orient_flag == 0:
|
||||
fit_sugars(dummy_res_A,old)
|
||||
elif double_stranded_bool == False:
|
||||
fit_sugars(dummy_res_A, old)
|
||||
else:
|
||||
raise pymol.CmdException("Double stranded bool was not provided to move_new_res")
|
||||
|
||||
dummy_fragment_transform = _self.get_object_matrix(dummy_fragment)
|
||||
_self.transform_object(full_frag, dummy_fragment_transform)
|
||||
_self.delete(dummy_fragment)
|
||||
|
||||
def check_dummy_oriention(old, dummy_res_A, *, _self=cmd):
|
||||
dummy_orient = _prefix + "_dummy_orient"
|
||||
orient_flag = 0
|
||||
orient_flag = _self.select(dummy_orient, f"({old} & name O4') within 1.0 of ({dummy_res_A} & name O4')")
|
||||
|
||||
return orient_flag
|
||||
|
||||
class NascentNucAcidInfo:
|
||||
def __init__(self, fragment_name, nuc_acid, nuc_type, form, dbl_helix):
|
||||
self.fragment_name = fragment_name
|
||||
self.nuc_acid = nuc_acid
|
||||
self.nuc_type = nuc_type
|
||||
self.form = form
|
||||
self.dbl_helix = dbl_helix
|
||||
|
||||
def attach_O5_phosphate(_self=cmd):
|
||||
if "pk1" not in _self.get_names("selections"):
|
||||
raise pymol.CmdException("Selection must be pk1 to attach O5' phosphate")
|
||||
|
||||
print("This building selection has an unphosphorylated O5' end.")
|
||||
attach_fragment("pk1","phosphite",4,0)
|
||||
# Initailize selection strings
|
||||
P_center = _prefix + "_P_center"
|
||||
H_extra = _prefix + "_H_extra"
|
||||
O_one = _prefix + "_O_one"
|
||||
O_two = _prefix + "_O_two"
|
||||
O_three = _prefix + "_O_three"
|
||||
|
||||
# Selection
|
||||
_self.select(P_center,"n. P01")
|
||||
_self.select(H_extra, f"h. and bound_to {P_center} or n. H02")
|
||||
_self.select(O_one, f"n. O01 and bound_to {P_center}")
|
||||
_self.select(O_two, f"n. O02 and bound_to {P_center}")
|
||||
_self.select(O_three, f"n. O03 and bound_to {P_center}")
|
||||
|
||||
# Removing unnecessary atoms
|
||||
_self.remove(H_extra)
|
||||
_self.remove(O_one)
|
||||
|
||||
# Fix bonding
|
||||
_self.unbond(P_center,O_three)
|
||||
_self.bond(P_center,O_three,1)
|
||||
_self.unbond(P_center,O_two)
|
||||
_self.bond(P_center,O_two,2)
|
||||
|
||||
# Rename P correctly
|
||||
_self.alter(P_center,"name = 'P'")
|
||||
|
||||
# Set Pk1 correctly
|
||||
_self.select("pk1",P_center)
|
||||
|
||||
def check_DNA_base_pair(sele_oppo_atom, selection, *, _self=cmd):
|
||||
base_pair_dist = _prefix + "_base_pair_dist"
|
||||
base_pair_bool = 0
|
||||
tmp_last_resn = iterate_to_list(selection,"resn")
|
||||
tmp_last_resn_oppo = iterate_to_list(sele_oppo_atom,"resn")
|
||||
|
||||
if len(tmp_last_resn_oppo) != 0:
|
||||
last_resn = str(tmp_last_resn[0])
|
||||
last_resn_oppo = str(tmp_last_resn_oppo[0])
|
||||
|
||||
if (_oneNA_base_pair[last_resn] == last_resn_oppo and
|
||||
_self.select(base_pair_dist, f"(byres {sele_oppo_atom}) within 3.5 of (byres {selection})") != 0):
|
||||
base_pair_bool = 1
|
||||
else:
|
||||
base_pair_bool = 0
|
||||
else:
|
||||
print("check_DNA_base_pair has no opposing residue to check")
|
||||
return base_pair_bool
|
||||
|
||||
def get_chains_oppo (chain, tmp_connect, *, _self=cmd):
|
||||
models = iterate_to_list(tmp_connect,"model")
|
||||
|
||||
close_chains = []
|
||||
close_chains = _self.get_chains(f"({models[0]}) within 15.0 of {tmp_connect}")
|
||||
close_chains = [c for c in close_chains if c != chain]
|
||||
|
||||
return close_chains
|
||||
|
||||
def get_new_chain (chain, tmp_connect, *, _self=cmd):
|
||||
models = iterate_to_list(tmp_connect,"model")
|
||||
model_chains = _self.get_chains(models[0])
|
||||
search_chain_flag = 0
|
||||
|
||||
if len(model_chains) != 0:
|
||||
last_chain = f"{model_chains[-1]}"
|
||||
last_chain_front = last_chain[:-1]
|
||||
last_chain_back = last_chain[-1]
|
||||
if last_chain_back != 'Z' and last_chain_back != 'z':
|
||||
new_chain_back = chr(ord(last_chain_back)+1)
|
||||
elif last_chain_back == 'Z':
|
||||
new_chain_back = "ZA"
|
||||
print("Z chain was detected. New chain will append A")
|
||||
else:
|
||||
new_chain_back = "za"
|
||||
print("z chain was detected. New chain will append a")
|
||||
else:
|
||||
new_chain_back = "A"
|
||||
|
||||
chain_oppo = last_chain_front + new_chain_back
|
||||
return chain_oppo
|
||||
|
||||
def check_valid_attachment(nascent, atom_selection_name, selection, resv, *, _self=cmd):
|
||||
atom_selection_name_partner = "O3'" if atom_selection_name == "P" else "P"
|
||||
atom_sele = _prefix + "atom_sele"
|
||||
_self.select(atom_sele, f"{selection}")
|
||||
bound = _prefix + "atom_sele_bound"
|
||||
if _self.count_atoms(f"(bound_to {atom_sele}) & name {atom_selection_name_partner}") != 0:
|
||||
_self.delete(tmp_wild)
|
||||
raise pymol.CmdException(f"{atom_selection_name} already bonded!")
|
||||
|
||||
def bond_single_stranded(tmp_editor, object, chain, resv, last_resi_sele, atom_selection_name, atom_name_oppo, *, _self=cmd):
|
||||
"""
|
||||
Forms a bond between the last atoms on selected structure and the newly created fragment
|
||||
:param tmp_editor: (str) Object representing the newly created fragment
|
||||
:param object: (str) object/model name of the selected structure
|
||||
:param chain: (str) Chain ID
|
||||
:param resv: (int)
|
||||
:param last_resi_sele: (str) The selection string of the selected residue
|
||||
:param atom_selection_name: (str) Name of the atom selected
|
||||
:param atom_name_oppo: (str) Name of the corresponding opposing atom
|
||||
"""
|
||||
object_fuse = _prefix + f"_{chain}_fuse"
|
||||
object_connect = _prefix + f"_{chain}_con"
|
||||
|
||||
print("The program did not detect a double stranded structure, so the opposing residue will not be attached.")
|
||||
|
||||
# Select and fuse
|
||||
_self.select(object_fuse, f"{last_resi_sele} & name {atom_selection_name}")
|
||||
_self.fuse(f"{tmp_editor} & chain {chain} & name {atom_name_oppo}", object_fuse, mode=3)
|
||||
_self.select(object_connect, f"{last_resi_sele} & name {atom_selection_name} & chain {chain}")
|
||||
|
||||
# Target is on the new fragment
|
||||
object_bond_target = _prefix + f"_{chain}_con_target"
|
||||
if (_self.select(object_bond_target, f"{object} & resi \\{resv} & name {atom_name_oppo} & chain {chain}") == 1):
|
||||
bond_dist = _prefix + "_bond_dist"
|
||||
if (_self.select(bond_dist, f"{object_connect} within 3.0 of {object_bond_target}") != 0):
|
||||
_self.bond(object_connect, object_bond_target)
|
||||
else:
|
||||
print("Identified bond targets were too far apart, so this will not be bound")
|
||||
else:
|
||||
print("More than one bond target was identified, so this will not be bound")
|
||||
|
||||
def bond_double_stranded(tmp_editor, object, chain, chain_oppo, resv, resv_oppo, last_resi_sele, prev_oppo_res, atom_selection_name,
|
||||
atom_name_oppo, *, _self=cmd):
|
||||
"""
|
||||
Forms a bond between the last atoms on selected structure and the newly created fragment
|
||||
:param tmp_editor: (str) Object representing the newly created fragment
|
||||
:param object: (str) object/model name of the selected structure
|
||||
:param chain: (str) Chain ID
|
||||
:param chain_oppo: (str) Opposing chain ID
|
||||
:param resv: (int)
|
||||
:param resv_oppo: (int)
|
||||
:param last_resi_sele: (str) The selection string of the selected residue
|
||||
:param prev_oppo_res: (str) The selection string of the previous opposing residue
|
||||
:param atom_selection_name: (str) Name of the atom selected
|
||||
:param atom_name_oppo: (str) Name of the corresponding opposing atom
|
||||
"""
|
||||
object_fuse = _prefix + f"_{chain}_fuse"
|
||||
object_oppo_fuse = _prefix + f"_oppo_{chain_oppo}_fuse"
|
||||
object_connect = _prefix + f"_{chain}_con"
|
||||
object_oppo_connect = _prefix + f"_oppo_{chain_oppo}_con"
|
||||
|
||||
# Target is on the new fragment
|
||||
object_bond_target = _prefix + f"_{chain}_con_target"
|
||||
object_oppo_bond_target = _prefix + f"_oppo_{chain_oppo}_con_target"
|
||||
|
||||
# Select and fuse
|
||||
_self.select(object_fuse, f"{last_resi_sele} & name {atom_selection_name}")
|
||||
_self.select(object_oppo_fuse, f"{prev_oppo_res} & name {atom_name_oppo}")
|
||||
_self.fuse(f"{tmp_editor} & chain {chain} & name {atom_name_oppo}", object_fuse, mode=3)
|
||||
|
||||
if ((_self.select(object_bond_target, f"{object} & resi \\{resv} & name {atom_name_oppo} & chain {chain}") == 1) and
|
||||
(_self.select(object_connect, f"{last_resi_sele} & name {atom_selection_name} & chain {chain}") == 1)):
|
||||
bond_dist = _prefix + "_bond_dist"
|
||||
if (_self.select(bond_dist, f"{object_connect} within 3.0 of {object_bond_target}") != 0):
|
||||
_self.bond(object_connect, object_bond_target)
|
||||
else:
|
||||
print("Identified bond targets were too far apart, so this will not be bound")
|
||||
else:
|
||||
print("More than one bond target was found on selected chain, so this will not be bound.")
|
||||
|
||||
if ((_self.select(object_oppo_bond_target, f"{object} & resi \\{resv_oppo} & name {atom_selection_name} & chain {chain_oppo}") == 1) and
|
||||
(_self.select(object_oppo_connect, f"{prev_oppo_res} & name {atom_name_oppo} & chain {chain_oppo}") == 1)):
|
||||
bond_dist = _prefix + "_bond_dist"
|
||||
if (_self.select(bond_dist, f"{object_oppo_connect} within 3.0 of {object_oppo_bond_target}") != 0):
|
||||
_self.bond(object_oppo_connect, object_oppo_bond_target)
|
||||
else:
|
||||
print("Identified bond targets were too far apart, so this will not be bound")
|
||||
else:
|
||||
print("More than one bond target was found on opposing chain, so this will not be bound.")
|
||||
|
||||
def attach_nuc_acid(selection, nuc_acid, nuc_type, object= "", form ="B",
|
||||
dbl_helix=True, *, _self=cmd):
|
||||
"""
|
||||
Creates new nuc acid attached to existing PDB structure
|
||||
:param selection: (str) selection of picked nascent chain (or nothing)
|
||||
:param nuc_acid: (str) appending nucleic acid
|
||||
:param nuc_type: (str) sugar type of nucleic acid
|
||||
:param object: (str) name of appending nucleobase
|
||||
:param form: (str) DNA structure form: A, B, or Z
|
||||
:param dbl_helix: (bool) flag for double-strandedness
|
||||
"""
|
||||
original_sele = _prefix + "_original_sele"
|
||||
_self.select(original_sele,selection)
|
||||
|
||||
if nuc_type == "RNA" and form != 'A':
|
||||
form = 'A'
|
||||
dbl_helix = False
|
||||
|
||||
nascent = NascentNucAcidInfo(nuc_acid + form, nuc_acid, nuc_type, form, dbl_helix)
|
||||
nuc_acid_partner_temp = _base_pair[nuc_type][nuc_acid].lower()
|
||||
nascent_partner = NascentNucAcidInfo(nuc_acid_partner_temp + form, nuc_acid_partner_temp,
|
||||
nuc_type, form, dbl_helix)
|
||||
|
||||
if _self.cmd.count_atoms(selection) == 0:
|
||||
if object == "":
|
||||
object = nuc_acid
|
||||
|
||||
if dbl_helix:
|
||||
frag_string = nascent.nuc_acid + "_"+ _base_pair["DNA"][nascent.nuc_acid] + nascent.form
|
||||
_self.fragment(frag_string,object)
|
||||
elif not dbl_helix:
|
||||
_self.fragment(nascent.fragment_name, object, origin=0)
|
||||
_self.alter(object, f"segi='A';chain='A';resv=1")
|
||||
rename_three_to_one(nascent.nuc_acid, object, nascent.nuc_type)
|
||||
|
||||
if nascent.nuc_type == "RNA":
|
||||
add2pO(object, nascent.nuc_acid, 1)
|
||||
|
||||
if _self.count_atoms(f"{object} & segi A & name P"):
|
||||
_self.edit(f"{object} & segi A & name P")
|
||||
elif _self.count_atoms(f"{object} & segi A & name O3'"):
|
||||
_self.edit(f"{object} & segi A & name O3'")
|
||||
_self.edit(f"{object} & segi A & name O3'")
|
||||
_self.select("pk1", f"{object} & name O3' & chain A")
|
||||
elif _self.select(tmp_connect, f'{selection}') == 1:
|
||||
chain, name, object = iterate_to_list(selection,"chain,name,model")[0]
|
||||
|
||||
if name == "O5'":
|
||||
attach_O5_phosphate()
|
||||
name = "P"
|
||||
# FIXME Don't use `selection` as the name, or ensure that it's
|
||||
# indeed just a name and not a selection expression.
|
||||
_self.select(selection,f"name P & bound_to {original_sele}")
|
||||
_self.delete("_pkdihe")
|
||||
_self.select(tmp_connect,f"{selection}")
|
||||
|
||||
if name == "P" or name == "O3'":
|
||||
extend_nuc_acid(nascent, nascent_partner, selection, object, name, chain, _self=_self)
|
||||
|
||||
_self.select("pk1","_tmp_editor_new_selection")
|
||||
else:
|
||||
_self.delete(tmp_wild)
|
||||
raise pymol.CmdException("invalid connection point: must be one atom, name O3' or P")
|
||||
|
||||
_self.show("cartoon", f"byobject {selection}")
|
||||
_self.delete(tmp_wild)
|
||||
|
||||
def extend_nuc_acid(nascent, nascent_partner, selection,
|
||||
object, atom_selection_name, chain = 'A', *, _self=cmd):
|
||||
"""
|
||||
Creates new nuc acid (or pair) or attaches to PDB chain
|
||||
:param nascent: (NascentNucAcidInfo) appending nucleic acid
|
||||
:param nascent_partner: (NascentNucAcidInfo) partner of appending nucleic acid
|
||||
:param selection: (str) selection string of selected residue
|
||||
:param object: (str) object of selected residue
|
||||
:param atom_selection_name: (str) O3' or P
|
||||
:param chain: (str) chain ID
|
||||
|
||||
FIXME Eliminate `tmp_connect` pre-condition (or at least document it
|
||||
properly). Looks like `tmp_connect` must be the same as `selection`?
|
||||
"""
|
||||
# Making a temporary selection
|
||||
original_sele = _prefix + "_original_sele"
|
||||
_self.select(original_sele,selection)
|
||||
|
||||
# Alter the segi to match the chain selection
|
||||
chain_sele = "_chain_sele"
|
||||
_self.select(chain_sele, f"chain {chain}")
|
||||
_self.alter(chain_sele,f"segi = '{chain}'")
|
||||
_self.delete(chain_sele)
|
||||
|
||||
if not nascent.dbl_helix:
|
||||
frag_string = nascent.fragment_name
|
||||
_self.fragment(frag_string, tmp_editor, origin=0)
|
||||
rename_three_to_one(nascent.nuc_acid, tmp_editor, nascent.nuc_type)
|
||||
elif nascent.dbl_helix:
|
||||
frag_string = nascent.nuc_acid + "_"+ _base_pair["DNA"][nascent.nuc_acid] + nascent.form
|
||||
_self.fragment(frag_string, tmp_editor, origin=0)
|
||||
else:
|
||||
raise pymol.CmdException("No helix state selected")
|
||||
|
||||
if _self.count_atoms(f"name {atom_selection_name}", domain=tmp_connect):
|
||||
tmp_resv = iterate_to_list(tmp_connect,"resv")
|
||||
|
||||
# Assign last res before adjustment is made
|
||||
last_resv = int(tmp_resv[0])
|
||||
|
||||
if atom_selection_name == "O3'":
|
||||
tmp_resv[0] = str(tmp_resv[0] + 1)
|
||||
elif atom_selection_name == "P":
|
||||
tmp_resv[0] = str(tmp_resv[0] - 1)
|
||||
else:
|
||||
raise pymol.CmdException("Something went wrong with resv loop in extend_nuc")
|
||||
|
||||
# Resv and resi assignment and testing
|
||||
resv = int(tmp_resv[0])
|
||||
resi = resv
|
||||
|
||||
check_valid_attachment(nascent, atom_selection_name, selection, resv)
|
||||
|
||||
reverse = False if atom_selection_name == "O3'" else True
|
||||
|
||||
last_resi_sele = _prefix + "_last_resi"
|
||||
_self.select(last_resi_sele, f"(byobject {selection}) & chain {chain} & resi \\{last_resv}")
|
||||
|
||||
if not nascent.dbl_helix:
|
||||
_self.alter(tmp_editor, f"chain='{chain}';segi='{chain}';resi=tmp_resv[0]",
|
||||
space={'tmp_resv': tmp_resv})
|
||||
# Set parameters so that move_new_res can be used
|
||||
prev_oppo_res = "none"
|
||||
double_stranded_bool = False
|
||||
|
||||
# Move and fuse
|
||||
move_new_res(frag_string, tmp_editor, last_resi_sele, prev_oppo_res, double_stranded_bool, nascent.form, antisense=reverse)
|
||||
_self.fuse(f"{tmp_editor} & name P", tmp_connect, mode=3)
|
||||
|
||||
_self.select(tmp_domain, "byresi (pk1 | pk2)")
|
||||
|
||||
if atom_selection_name == "O3'":
|
||||
_self.bond(f"{tmp_domain} & resi \\{last_resv} & name O3'",
|
||||
f"{tmp_domain} & resi \\{resv} & name P")
|
||||
elif atom_selection_name == "P":
|
||||
_self.bond(f"{tmp_domain} & resi \\{resv} & name O3'",
|
||||
f"{tmp_domain} & resi \\{last_resv} & name P")
|
||||
|
||||
elif nascent.dbl_helix:
|
||||
#Initialize strings for selection
|
||||
prev_oppo_res = _prefix + "_prev_oppo_res"
|
||||
end_oppo_atom = _prefix + "_end_oppo_atom"
|
||||
|
||||
#Initialize a double stranded bool for the detection of existing double strand
|
||||
double_stranded_bool = False
|
||||
base_pair_result = 0
|
||||
|
||||
# Get oppo information
|
||||
atom_name_oppo = "O3'" if atom_selection_name == "P" else "P"
|
||||
|
||||
# Get opposing chains and iterate looking for basepairs
|
||||
chains_oppo = get_chains_oppo(chain, tmp_connect)
|
||||
for tmp_chain_oppo in chains_oppo:
|
||||
# Create selection strings
|
||||
last_oppo_atom = _prefix + "_last_oppo_atom"
|
||||
first_oppo_atom = _prefix + "_first_oppo_atom"
|
||||
same_oppo_res = _prefix + "_same_oppo_atom"
|
||||
|
||||
_self.select(last_oppo_atom, f"last (({object} and chain {tmp_chain_oppo} and (not {tmp_editor})) and polymer.nuc)")
|
||||
_self.select(first_oppo_atom, f"first (({object} and chain {tmp_chain_oppo} and (not {tmp_editor})) and polymer.nuc)")
|
||||
base_pair_last = check_DNA_base_pair(last_oppo_atom, original_sele)
|
||||
base_pair_first = check_DNA_base_pair(first_oppo_atom, original_sele)
|
||||
|
||||
same_base_flag = _self.select(same_oppo_res, f"({object} and (byres {last_oppo_atom}) and (byres {first_oppo_atom}))")
|
||||
|
||||
if ((base_pair_first == 1 or base_pair_last == 1) and base_pair_result == 1 and same_base_flag == 0):
|
||||
print("Multiple residues meet base pairing requirements. Building as if no opposing strand detected.")
|
||||
base_pair_result = 0
|
||||
break
|
||||
elif (base_pair_first == 1 and base_pair_last ==1):
|
||||
chain_oppo = tmp_chain_oppo
|
||||
base_pair_result = 1
|
||||
|
||||
check_DNA_base_pair(first_oppo_atom, original_sele)
|
||||
atoms_first = iterate_to_list(_prefix + "_base_pair_dist","name")
|
||||
check_DNA_base_pair(last_oppo_atom, original_sele)
|
||||
atoms_last = iterate_to_list(_prefix + "_base_pair_dist", "name")
|
||||
|
||||
if (len(atoms_first) > len(atoms_last)):
|
||||
_self.select(end_oppo_atom,first_oppo_atom)
|
||||
else:
|
||||
_self.select(end_oppo_atom,last_oppo_atom)
|
||||
|
||||
elif (base_pair_last == 1):
|
||||
base_pair_result = 1
|
||||
chain_oppo = tmp_chain_oppo
|
||||
_self.select(end_oppo_atom,last_oppo_atom)
|
||||
elif (base_pair_first == 1):
|
||||
base_pair_result = 1
|
||||
chain_oppo = tmp_chain_oppo
|
||||
_self.select(end_oppo_atom,first_oppo_atom)
|
||||
else:
|
||||
print("No based pair was found on chain ", tmp_chain_oppo)
|
||||
|
||||
if (base_pair_result == 1):
|
||||
double_stranded_bool = True
|
||||
# Confirm base pair result and set _base_pair_dist
|
||||
check_DNA_base_pair(end_oppo_atom, original_sele)
|
||||
|
||||
tmp_last_resv_oppo = iterate_to_list(end_oppo_atom,"resv")
|
||||
if len(tmp_last_resv_oppo) == 1:
|
||||
last_resv_oppo = tmp_last_resv_oppo[0]
|
||||
else:
|
||||
# This is arbitrary and may be changed. Using neg selected resv for now.
|
||||
last_resv_oppo = -last_resv
|
||||
|
||||
_self.select(prev_oppo_res, "byres (_tmp_editor_base_pair_dist)")
|
||||
elif (base_pair_result == 0):
|
||||
last_resv_oppo = -last_resv
|
||||
chain_oppo = get_new_chain(chain, tmp_connect)
|
||||
else:
|
||||
raise pymol.CmdException("Base pairing result is not returning 0 or 1")
|
||||
|
||||
# Alter the opposing segi to match chain
|
||||
chain_oppo_sele = "_chain_oppo_sele"
|
||||
_self.select(chain_oppo_sele, f"chain {chain_oppo}")
|
||||
_self.alter(chain_oppo_sele,f"segi = '{chain_oppo}'")
|
||||
_self.delete(chain_oppo_sele)
|
||||
|
||||
# Get oppo resv value based on selection
|
||||
if atom_selection_name == "O3'":
|
||||
resv_oppo = last_resv_oppo - 1
|
||||
elif atom_selection_name == "P":
|
||||
resv_oppo = last_resv_oppo + 1
|
||||
|
||||
# If picking O3', check O5' phosphate based on info found in this loop
|
||||
if (atom_selection_name == "O3'" and double_stranded_bool == True):
|
||||
tmp_phosphate_check = _prefix + "_phosphate_check"
|
||||
|
||||
if(_self.select(tmp_phosphate_check, f"{prev_oppo_res} and name P ") == 0):
|
||||
_self.select("pk1", f"{prev_oppo_res} and name O5'")
|
||||
_self.select("pk2", f"{prev_oppo_res} and name O5'") # This is so attach frag has a pk2
|
||||
|
||||
attach_O5_phosphate()
|
||||
|
||||
_self.unpick()
|
||||
_self.select("pk1",original_sele)
|
||||
_self.select(tmp_connect,f"{selection}")
|
||||
_self.select(prev_oppo_res, f"byres ({prev_oppo_res})")
|
||||
|
||||
if (_self.select(tmp_phosphate_check, f"{prev_oppo_res} and name P")==1):
|
||||
print("Phosphate has been successfully added")
|
||||
|
||||
# Move the fragment
|
||||
move_new_res(frag_string, tmp_editor, last_resi_sele, prev_oppo_res, double_stranded_bool, nascent.form, antisense=reverse)
|
||||
|
||||
# Alter residues created
|
||||
_self.alter("_tmp_editor_newresA",f"chain='{chain}';segi='{chain}';resv={resv}")
|
||||
_self.alter("_tmp_editor_newresB",f"chain='{chain_oppo}';segi='{chain_oppo}';resv={resv_oppo}")
|
||||
|
||||
# Connect the created residues
|
||||
if double_stranded_bool == True:
|
||||
bond_double_stranded(tmp_editor, object, chain, chain_oppo, resv, resv_oppo, last_resi_sele, prev_oppo_res,
|
||||
atom_selection_name, atom_name_oppo)
|
||||
elif double_stranded_bool == False:
|
||||
bond_single_stranded(tmp_editor, object, chain, resv, last_resi_sele, atom_selection_name, atom_name_oppo)
|
||||
else:
|
||||
raise pymol.CmdException("double_stranded_bool is not returning True or False")
|
||||
|
||||
if nascent.nuc_type == "RNA":
|
||||
add2pO(tmp_domain, nascent.nuc_acid, resv)
|
||||
|
||||
new_term = _prefix + "_new_selection"
|
||||
_self.select(new_term, f"{object} & chain {chain} & resi \\{resv} & name {atom_selection_name}")
|
||||
#_self.edit(new_term)
|
||||
|
||||
def fab(input,name=None,mode='peptide',resi=1,chain='',segi='',state=-1,
|
||||
dir=1,hydro=-1,ss=0,async_=0,quiet=1,_self=cmd, **kwargs):
|
||||
'''
|
||||
@@ -422,6 +1077,73 @@ EXAMPLE
|
||||
r = DEFAULT_SUCCESS
|
||||
return r
|
||||
|
||||
def fnab(input, name=None, mode="DNA", form="B", dbl_helix=1, *, _self=cmd):
|
||||
"""
|
||||
DESCRIPTION
|
||||
|
||||
Builds a nucleotide acid from sequence
|
||||
|
||||
Fragments provided by:
|
||||
Lu, Xiang-Jun, Olson, Wilma K. 3DNA: a software package for the analysis,
|
||||
rebuilding and visualization of three-dimensional nucleic acid structures,
|
||||
Nucleic Acids Research, 32, W667-W675 (2004).
|
||||
|
||||
USAGE
|
||||
|
||||
fnab input [, name [, type [, form [, dbl_helix ]]]]
|
||||
|
||||
ARGUMENTS
|
||||
|
||||
input = str: Sequence as an array of one letter codes
|
||||
|
||||
name = str: Name of the object to create {default: obj}
|
||||
|
||||
mode = str: "DNA" or "RNA"
|
||||
|
||||
form = str: "A" or "B"
|
||||
|
||||
dbl_helix = bool (0/1): flag for using double helix in DNA
|
||||
|
||||
EXAMPLE
|
||||
|
||||
fnab ATGCGATAC
|
||||
fnab ATGCGATAC, name=myDNA, mode=DNA, form=B, dbl_helix=1
|
||||
fnab AAUUUUCCG, mode=RNA
|
||||
"""
|
||||
_self.unpick()
|
||||
|
||||
if name is None:
|
||||
name = _self.get_unused_name(prefix="obj")
|
||||
elif name in _self.get_names('all'):
|
||||
name = _self.get_unused_name(prefix=name)
|
||||
|
||||
dbl_helix = int(dbl_helix) > 0
|
||||
|
||||
mode = mode.upper()
|
||||
if mode not in ('DNA', 'RNA'):
|
||||
raise pymol.CmdException("\"mode\" must be \"DNA\" or \"RNA\" only.")
|
||||
if mode == "RNA" and dbl_helix != 0:
|
||||
print ("Double helix RNA building is not currently supported.")
|
||||
dbl_helix = 0
|
||||
|
||||
#first pass for error checking
|
||||
for oneNA in input:
|
||||
oneNA = oneNA.upper()
|
||||
threeNA = _oneNA_to_threeNA.get(oneNA)
|
||||
if threeNA is None:
|
||||
raise pymol.CmdException("\"%s\" not of %s type..." % (oneNA, mode))
|
||||
|
||||
if threeNA not in _base_pair[mode]:
|
||||
raise pymol.CmdException("\"%s\" not of %s type..." % (oneNA, mode))
|
||||
|
||||
for oneNA in input:
|
||||
oneNA = oneNA.upper()
|
||||
threeNA = _oneNA_to_threeNA[oneNA]
|
||||
attach_nuc_acid(selection="?pk1", nuc_acid=threeNA, nuc_type=mode,
|
||||
object=name, form=form, dbl_helix=dbl_helix)
|
||||
_self.unpick()
|
||||
return DEFAULT_SUCCESS
|
||||
|
||||
def build_peptide(sequence,_self=cmd): # legacy
|
||||
for aa in sequence:
|
||||
attach_amino_acid("pk1",_aa_codes[aa])
|
||||
|
||||
@@ -83,6 +83,7 @@ def get_command_keywords(self_cmd=cmd):
|
||||
'extract' : [ self_cmd.extract , 0 , 0 , '' , parsing.STRICT ],
|
||||
'exec' : [ self_cmd.python_help , 0 , 0 , '' , parsing.PYTHON ],
|
||||
'fab' : [ self_cmd.fab , 0 , 0 , '' , parsing.STRICT ],
|
||||
'fnab' : [ self_cmd.fnab , 0 , 0 , '' , parsing.STRICT ],
|
||||
'feedback' : [ self_cmd.feedback , 0, 0 , '' , parsing.STRICT ],
|
||||
'fetch' : [ self_cmd.fetch , 0, 0 , '' , parsing.STRICT ],
|
||||
'fit' : [ self_cmd.fit , 0 , 0 , '' , parsing.STRICT ],
|
||||
|
||||
@@ -3,7 +3,6 @@ Testing PyQt Fab & Builder
|
||||
"""
|
||||
from pymol import cmd, testing
|
||||
|
||||
@testing.requires('incentive')
|
||||
@testing.requires_version('2.3')
|
||||
class TestNucBuilder(testing.PyMOLTestCase):
|
||||
def testCanInit(self):
|
||||
|
||||
Reference in New Issue
Block a user