1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53
|
from ost import io
from promod3 import loop, scoring, modelling
# let's load a crambin structure from the pdb
crambin = io.LoadPDB('data/1CRN.pdb')
SEQRES = ''.join([r.one_letter_code for r in crambin.residues])
# this is the sequence we want to remodel
loop_seq = SEQRES[23:31]
# let's define the stem residues
n_stem = crambin.residues[23]
c_stem = crambin.residues[30]
# we use the StructureDB as source for structural information
structure_db = loop.LoadStructureDB()
# the FragDB allows to access the StructureDB based on geometric
# features of the loop stem residue
frag_db = loop.LoadFragDB()
# the LoopCandidates allow to handle several loops at once
# we directly want to find potential loop candidates from the
# previously loaded databases
loop_candidates = modelling.LoopCandidates.FillFromDatabase(\
n_stem, c_stem, loop_seq, frag_db, structure_db)
# candidates usually don't match exactly the required stem coords.
# CCD (Cyclic Coordinate Descent) is one way to enforce this match.
loop_candidates.ApplyCCD(n_stem, c_stem)
# setup backbone scorer with clash and cbeta scoring
score_env = scoring.BackboneScoreEnv(SEQRES)
score_env.SetInitialEnvironment(crambin)
scorer = scoring.BackboneOverallScorer()
scorer["cbeta"] = scoring.LoadCBetaScorer()
scorer["clash"] = scoring.ClashScorer()
scorer.AttachEnvironment(score_env)
# the scorer can then be used on the LoopCandidates object to
# calculate desired scores (here: cbeta & clash, start resnum = 24)
bb_scores = modelling.ScoreContainer()
loop_candidates.CalculateBackboneScores(bb_scores, scorer, score_env,
["cbeta", "clash"], 24)
# we finally perform a weighted linear combination of the scores,
# pick the best one, insert it into our structure and save it
weights = {"cbeta": 1, "clash": 1}
scores = bb_scores.LinearCombine(weights)
min_score = min(scores)
min_candidate = scores.index(min_score)
bb_list = loop_candidates[min_candidate]
bb_list.InsertInto(crambin.chains[0], n_stem.GetNumber())
io.SavePDB(crambin, "modified_crambin.pdb")
|