3 """@namespace IMP.pmi1.tools
4 Miscellaneous utilities.
7 from __future__
import print_function
15 from math
import log,pi,sqrt,exp
22 from collections
import defaultdict
24 from collections
import OrderedDict
26 from IMP.pmi1._compat_collections
import OrderedDict
28 def _add_pmi_provenance(p):
29 """Tag the given particle as being created by the current version of PMI."""
33 location=
"https://integrativemodeling.org")
36 def _get_restraint_set_keys():
37 if not hasattr(_get_restraint_set_keys,
'pmi_rs_key'):
38 _get_restraint_set_keys.pmi_rs_key =
IMP.ModelKey(
"PMI restraints")
39 _get_restraint_set_keys.rmf_rs_key =
IMP.ModelKey(
"RMF restraints")
40 return (_get_restraint_set_keys.pmi_rs_key,
41 _get_restraint_set_keys.rmf_rs_key)
43 def _add_restraint_sets(model, mk, mk_rmf):
46 model.add_data(mk, rs)
47 model.add_data(mk_rmf, rs_rmf)
51 """Add a PMI restraint to the model.
52 Since Model.add_restraint() no longer exists (in modern IMP restraints
53 should be added to a ScoringFunction instead) store them instead in
54 a RestraintSet, and keep a reference to it in the Model.
56 If `add_to_rmf` is True, also add the restraint to a separate list
57 of restraints that will be written out to RMF files (by default, most
58 PMI restraints are not)."""
59 mk, mk_rmf = _get_restraint_set_keys()
60 if model.get_has_data(mk):
61 rs = IMP.RestraintSet.get_from(model.get_data(mk))
62 rs_rmf = IMP.RestraintSet.get_from(model.get_data(mk_rmf))
64 rs, rs_rmf = _add_restraint_sets(model, mk, mk_rmf)
65 rs.add_restraint(restraint)
67 rs_rmf.add_restraint(restraint)
70 """Get a RestraintSet containing all PMI restraints added to the model.
71 If `rmf` is True, return only the subset of these restraints that
72 should be written out to RMF files."""
73 mk, mk_rmf = _get_restraint_set_keys()
74 if not model.get_has_data(mk):
75 print(
"WARNING: no restraints added to model yet")
76 _add_restraint_sets(model, mk, mk_rmf)
78 return IMP.RestraintSet.get_from(model.get_data(mk_rmf))
80 return IMP.RestraintSet.get_from(model.get_data(mk))
83 """Collect timing information.
84 Add an instance of this class to outputobjects to get timing information
89 @param isdelta if True (the default) then report the time since the
90 last use of this class; if False, report cumulative time."""
91 self.starttime = time.clock()
93 self.isdelta = isdelta
95 def set_label(self, labelstr):
101 newtime = time.clock()
102 output[
"Stopwatch_" + self.label +
"_delta_seconds"] \
103 = str(newtime - self.starttime)
104 self.starttime = newtime
106 output[
"Stopwatch_" + self.label +
"_elapsed_seconds"] \
107 = str(time.clock() - self.starttime)
111 class SetupNuisance(object):
113 def __init__(self, m, initialvalue, minvalue, maxvalue, isoptimized=True):
117 nuisance.set_lower(minvalue)
119 nuisance.set_upper(maxvalue)
122 nuisance.set_is_optimized(nuisance.get_nuisance_key(), isoptimized)
123 self.nuisance = nuisance
125 def get_particle(self):
129 class SetupWeight(object):
131 def __init__(self, m, isoptimized=True):
134 self.weight.set_weights_are_optimized(isoptimized)
136 def get_particle(self):
140 class SetupSurface(object):
142 def __init__(self, m, center, normal, isoptimized=True):
145 self.surface.set_coordinates_are_optimized(isoptimized)
146 self.surface.set_normal_is_optimized(isoptimized)
148 def get_particle(self):
153 "If you use this class please let the PMI developers know.")
154 class ParticleToSampleFilter(object):
155 def __init__(self, sampled_objects):
156 self.sampled_objects=sampled_objects
159 def add_filter(self,filter_string):
160 self.filters.append(filter_string)
162 def get_particles_to_sample(self):
163 particles_to_sample={}
164 for so
in self.sampled_objects:
165 ps_dict=so.get_particles_to_sample()
167 for f
in self.filters:
169 if key
not in particles_to_sample:
170 particles_to_sample[key]=ps_dict[key]
172 particles_to_sample[key]+=ps_dict[key]
173 return particles_to_sample
175 class ParticleToSampleList(object):
179 self.dictionary_particle_type = {}
180 self.dictionary_particle_transformation = {}
181 self.dictionary_particle_name = {}
188 particle_transformation,
190 if not particle_type
in [
"Rigid_Bodies",
"Floppy_Bodies",
"Nuisances",
"X_coord",
"Weights",
"Surfaces"]:
191 raise TypeError(
"not the right particle type")
193 self.dictionary_particle_type[particle] = particle_type
194 if particle_type ==
"Rigid_Bodies":
195 if type(particle_transformation) == tuple
and len(particle_transformation) == 2
and type(particle_transformation[0]) == float
and type(particle_transformation[1]) == float:
196 self.dictionary_particle_transformation[
197 particle] = particle_transformation
198 self.dictionary_particle_name[particle] = name
200 raise TypeError(
"ParticleToSampleList: not the right transformation format for Rigid_Bodies, should be a tuple of floats")
201 elif particle_type ==
"Surfaces":
202 if type(particle_transformation) == tuple
and len(particle_transformation) == 3
and all(isinstance(x, float)
for x
in particle_transformation):
203 self.dictionary_particle_transformation[
204 particle] = particle_transformation
205 self.dictionary_particle_name[particle] = name
207 raise TypeError(
"ParticleToSampleList: not the right transformation format for Surfaces, should be a tuple of floats")
209 if type(particle_transformation) == float:
210 self.dictionary_particle_transformation[
211 particle] = particle_transformation
212 self.dictionary_particle_name[particle] = name
214 raise TypeError(
"ParticleToSampleList: not the right transformation format, should be a float")
216 def get_particles_to_sample(self):
218 for particle
in self.dictionary_particle_type:
219 key = self.dictionary_particle_type[
220 particle] +
"ParticleToSampleList_" + self.dictionary_particle_name[particle] +
"_" + self.label
223 self.dictionary_particle_transformation[particle])
228 class Variance(object):
230 def __init__(self, model, tau, niter, prot, th_profile, write_data=False):
233 self.write_data = write_data
237 particles = IMP.atom.get_by_type(prot, IMP.atom.ATOM_TYPE)
238 self.particles = particles
240 self.refpos = [
IMP.core.XYZ(p).get_coordinates()
for p
in particles]
241 self.model_profile = th_profile
243 def perturb_particles(self, perturb=True):
244 for i, p
in enumerate(self.particles):
245 newpos = array(self.refpos[i])
247 newpos += random.normal(0, self.tau, 3)
251 def get_profile(self):
252 model_profile = self.model_profile
253 p = model_profile.calculate_profile(self.particles, IMP.saxs.CA_ATOMS)
254 return array([model_profile.get_intensity(i)
for i
in
255 range(model_profile.size())])
257 def init_variances(self):
259 N = self.model_profile.size()
260 a = self.profiles[0][:]
262 self.V = self.m * self.m.T
263 self.normm = linalg.norm(self.m)
264 self.normV = linalg.norm(self.V)
266 def update_variances(self):
267 a = matrix(self.profiles[-1])
268 n = float(len(self.profiles))
269 self.m = a.T / n + (n - 1) / n * self.m
270 self.V = a.T * a + self.V
271 self.oldnormm = self.normm
272 self.oldnormV = self.normV
273 self.normm = linalg.norm(self.m)
274 self.normV = linalg.norm(self.V)
275 self.diffm = (self.oldnormm - self.normm) / self.oldnormm
276 self.diffV = (self.oldnormV - self.normV) / self.oldnormV
278 def get_direct_stats(self, a):
283 for q, I
in enumerate(prof):
288 Sigma = (matrix(a - m))
289 Sigma = Sigma.T * Sigma / (nprof - 1)
290 mi = matrix(diag(1. / m))
291 Sigmarel = mi.T * Sigma * mi
292 return m, V, Sigma, Sigmarel
294 def store_data(self):
295 if not os.path.isdir(
'data'):
297 profiles = matrix(self.profiles)
298 self.directm, self.directV, self.Sigma, self.Sigmarel = \
299 self.get_direct_stats(array(profiles))
300 directV = self.directV
302 save(
'data/profiles', profiles)
304 fl = open(
'data/profiles.dat',
'w')
305 for i, l
in enumerate(array(profiles).T):
306 self.model_profile.get_q(i)
309 fl.write(
'%s ' % (k - self.directm[i]))
312 fl = open(
'data/profiles_rel.dat',
'w')
313 for i, l
in enumerate(array(profiles).T):
314 self.model_profile.get_q(i)
317 fl.write(
'%s ' % ((k - self.directm[i]) / self.directm[i]))
319 save(
'data/m', self.directm)
320 save(
'data/V', self.directV)
322 save(
'data/Sigma', Sigma)
324 fl = open(
'data/Sigma.dat',
'w')
325 model_profile = self.model_profile
326 for i
in range(model_profile.size()):
327 qi = model_profile.get_q(i)
328 for j
in range(model_profile.size()):
329 qj = model_profile.get_q(j)
330 vij = self.Sigma[i, j]
331 fl.write(
'%s %s %s\n' % (qi, qj, vij))
334 fl = open(
'data/eigenvals',
'w')
335 for i
in linalg.eigvalsh(Sigma):
337 Sigmarel = self.Sigmarel
338 save(
'data/Sigmarel', Sigmarel)
340 fl = open(
'data/Sigmarel.dat',
'w')
341 model_profile = self.model_profile
342 for i
in range(model_profile.size()):
343 qi = model_profile.get_q(i)
344 for j
in range(model_profile.size()):
345 qj = model_profile.get_q(j)
346 vij = self.Sigmarel[i, j]
347 fl.write(
'%s %s %s\n' % (qi, qj, vij))
350 fl = open(
'data/eigenvals_rel',
'w')
351 for i
in linalg.eigvalsh(Sigmarel):
354 fl = open(
'data/mean.dat',
'w')
355 for i
in range(len(self.directm)):
356 qi = self.model_profile.get_q(i)
358 fl.write(
'%s ' % self.directm[i])
359 fl.write(
'%s ' % sqrt(self.Sigma[i, i]))
362 def try_chol(self, jitter):
365 linalg.cholesky(Sigma + matrix(eye(len(Sigma))) * jitter)
366 except linalg.LinAlgError:
367 print(
"Decomposition failed with jitter =", jitter)
369 print(
"Successful decomposition with jitter =", jitter)
372 self.profiles = [self.get_profile()]
374 for n
in range(self.niter):
375 self.perturb_particles()
376 self.profiles.append(self.get_profile())
388 def get_cov(self, relative=True):
398 number_of_cross_links=10,
399 ambiguity_probability=0.1,
400 confidence_score_range=[0,100],
401 avoid_same_particles=
False):
402 '''Return a random cross-link dataset as a string.
403 Every line is a residue pair, together with UniqueIdentifier
406 residue_pairs=get_random_residue_pairs(representation, resolution, number_of_cross_links, avoid_same_particles=avoid_same_particles)
409 cmin=float(min(confidence_score_range))
410 cmax=float(max(confidence_score_range))
414 for (name1, r1, name2, r2)
in residue_pairs:
415 if random.random() > ambiguity_probability:
417 score=random.random()*(cmax-cmin)+cmin
418 dataset+=str(name1)+
" "+str(name2)+
" "+str(r1)+
" "+str(r2)+
" "+str(score)+
" "+str(unique_identifier)+
"\n"
425 def get_cross_link_data(directory, filename, dist, omega, sigma,
426 don=
None, doff=
None, prior=0, type_of_profile=
"gofr"):
428 (distmin, distmax, ndist) = dist
429 (omegamin, omegamax, nomega) = omega
430 (sigmamin, sigmamax, nsigma) = sigma
433 with open(filen)
as xlpot:
434 dictionary = ast.literal_eval(xlpot.readline())
436 xpot = dictionary[directory][filename][
"distance"]
437 pot = dictionary[directory][filename][type_of_profile]
439 dist_grid =
get_grid(distmin, distmax, ndist,
False)
440 omega_grid = get_log_grid(omegamin, omegamax, nomega)
441 sigma_grid = get_log_grid(sigmamin, sigmamax, nsigma)
443 if not don
is None and not doff
is None:
465 def get_cross_link_data_from_length(length, xxx_todo_changeme3, xxx_todo_changeme4, xxx_todo_changeme5):
466 (distmin, distmax, ndist) = xxx_todo_changeme3
467 (omegamin, omegamax, nomega) = xxx_todo_changeme4
468 (sigmamin, sigmamax, nsigma) = xxx_todo_changeme5
470 dist_grid =
get_grid(distmin, distmax, ndist,
False)
471 omega_grid = get_log_grid(omegamin, omegamax, nomega)
472 sigma_grid = get_log_grid(sigmamin, sigmamax, nsigma)
478 def get_grid(gmin, gmax, ngrid, boundaries):
480 dx = (gmax - gmin) / float(ngrid)
481 for i
in range(0, ngrid + 1):
482 if(
not boundaries
and i == 0):
484 if(
not boundaries
and i == ngrid):
486 grid.append(gmin + float(i) * dx)
492 def get_log_grid(gmin, gmax, ngrid):
494 for i
in range(0, ngrid + 1):
495 grid.append(gmin * exp(float(i) / ngrid * log(gmax / gmin)))
503 example '"{ID_Score}" > 28 AND "{Sample}" ==
504 "%10_1%" OR ":Sample}" == "%10_2%" OR ":Sample}"
505 == "%10_3%" OR ":Sample}" == "%8_1%" OR ":Sample}" == "%8_2%"'
508 import pyparsing
as pp
510 operator = pp.Regex(
">=|<=|!=|>|<|==|in").setName(
"operator")
511 value = pp.QuotedString(
513 r"[+-]?\d+(:?\.\d*)?(:?[eE][+-]?\d+)?")
514 identifier = pp.Word(pp.alphas, pp.alphanums +
"_")
515 comparison_term = identifier | value
516 condition = pp.Group(comparison_term + operator + comparison_term)
518 expr = pp.operatorPrecedence(condition, [
519 (
"OR", 2, pp.opAssoc.LEFT, ),
520 (
"AND", 2, pp.opAssoc.LEFT, ),
523 parsedstring = str(expr.parseString(inputstring)) \
529 .replace(
"{",
"float(entry['") \
530 .replace(
"}",
"'])") \
531 .replace(
":",
"str(entry['") \
532 .replace(
"}",
"'])") \
533 .replace(
"AND",
"and") \
538 def open_file_or_inline_text(filename):
540 fl = open(filename,
"r")
542 fl = filename.split(
"\n")
549 for i
in range(0, len(prot0) - 1):
550 for j
in range(i + 1, len(prot0)):
553 drmsd += (dist0 - dist1) ** 2
555 return math.sqrt(drmsd / npairs)
560 def get_ids_from_fasta_file(fastafile):
562 with open(fastafile)
as ff:
571 this function works with plain hierarchies, as read from the pdb,
572 no multi-scale hierarchies
579 atom_type=IMP.atom.AT_CA)
587 print(
"get_closest_residue_position: exiting while loop without result")
589 p = sel.get_selected_particles()
594 print(
"get_closest_residue_position: got NO residues for hierarchy %s and residue %i" % (hier, resindex))
595 raise Exception(
"get_closest_residue_position: got NO residues for hierarchy %s and residue %i" % (
598 raise ValueError(
"got multiple residues for hierarchy %s and residue %i; the list of particles is %s" % (hier, resindex, str([pp.get_name()
for pp
in p])))
602 Get the particle of the terminal residue at the GIVEN resolution
603 (NOTE: not the closest resolution!).
604 To get the terminal residue at the closest resolution use:
605 particles=IMP.pmi1.tools.select_by_tuple(representation,molecule_name)
606 particles[0] and particles[-1] will be the first and last particles
607 corresponding to the two termini.
608 It is needed for instance to determine the last residue of a pdb.
609 @param hier hierarchy containing the terminal residue
610 @param terminus either 'N' or 'C'
611 @param resolution resolution to use.
617 resolution=resolution,
624 if termresidue
is None:
625 termresidue = max(residues)
627 elif max(residues) >= termresidue:
628 termresidue = max(residues)
630 elif terminus ==
"N":
631 if termresidue
is None:
632 termresidue = min(residues)
634 elif min(residues) <= termresidue:
635 termresidue = min(residues)
638 raise ValueError(
"terminus argument should be either N or C")
643 """Get XYZ coordinates of the terminal residue at the GIVEN resolution"""
649 Return the residue index gaps and contiguous segments in the hierarchy.
651 @param hierarchy hierarchy to examine
652 @param start first residue index
653 @param end last residue index
655 @return A list of lists of the form
656 [[1,100,"cont"],[101,120,"gap"],[121,200,"cont"]]
659 for n, rindex
in enumerate(range(start, end + 1)):
661 atom_type=IMP.atom.AT_CA)
663 if len(sel.get_selected_particles()) == 0:
667 rindexcont = start - 1
668 if rindexgap == rindex - 1:
674 gaps.append([rindex, rindex,
"gap"])
680 rindexgap = start - 1
682 if rindexcont == rindex - 1:
689 gaps.append([rindex, rindex,
"cont"])
700 def set_map_element(self, xvalue, yvalue):
701 self.map[xvalue] = yvalue
703 def get_map_element(self, invalue):
704 if type(invalue) == float:
708 dist = (invalue - x) * (invalue - x)
717 return self.map[minx]
718 elif type(invalue) == str:
719 return self.map[invalue]
721 raise TypeError(
"wrong type for map")
727 selection_arguments=
None,
729 name_is_ambiguous=
False,
733 representation_type=
None):
735 this function uses representation=SimplifiedModel
736 it returns the corresponding selected particles
737 representation_type="Beads", "Res:X", "Densities", "Representation", "Molecule"
740 if resolution
is None:
742 resolution_particles =
None
743 hierarchies_particles =
None
744 names_particles =
None
745 residue_range_particles =
None
746 residue_particles =
None
747 representation_type_particles =
None
749 if not resolution
is None:
750 resolution_particles = []
751 hs = representation.get_hierarchies_at_given_resolution(resolution)
755 if not hierarchies
is None:
756 hierarchies_particles = []
757 for h
in hierarchies:
762 if name_is_ambiguous:
763 for namekey
in representation.hier_dict:
766 representation.hier_dict[namekey])
767 elif name
in representation.hier_dict:
770 print(
"select: component %s is not there" % name)
772 if not first_residue
is None and not last_residue
is None:
774 residue_indexes=range(first_residue, last_residue + 1))
776 for p
in sel.get_selected_particles()]
778 if not residue
is None:
781 for p
in sel.get_selected_particles()]
783 if not representation_type
is None:
784 representation_type_particles = []
785 if representation_type ==
"Molecule":
786 for name
in representation.hier_representation:
787 for repr_type
in representation.hier_representation[name]:
788 if repr_type ==
"Beads" or "Res:" in repr_type:
789 h = representation.hier_representation[name][repr_type]
792 elif representation_type ==
"PDB":
793 for name
in representation.hier_representation:
794 for repr_type
in representation.hier_representation[name]:
795 if repr_type ==
"Res:" in repr_type:
796 h = representation.hier_representation[name][repr_type]
800 for name
in representation.hier_representation:
801 h = representation.hier_representation[
806 selections = [hierarchies_particles, names_particles,
807 residue_range_particles, residue_particles, representation_type_particles]
809 if resolution
is None:
810 selected_particles = set(allparticles)
812 selected_particles = set(resolution_particles)
816 selected_particles = (set(s) & selected_particles)
818 return list(selected_particles)
825 name_is_ambiguous=
False):
826 if isinstance(tupleselection, tuple)
and len(tupleselection) == 3:
828 name=tupleselection[2],
829 first_residue=tupleselection[0],
830 last_residue=tupleselection[1],
831 name_is_ambiguous=name_is_ambiguous)
832 elif isinstance(tupleselection, str):
835 name_is_ambiguous=name_is_ambiguous)
837 raise ValueError(
'you passed something bad to select_by_tuple()')
839 particles = IMP.pmi1.tools.sort_by_residues(particles)
844 """New tuple format: molname OR (start,stop,molname,copynum,statenum)
845 Copy and state are optional. Can also use 'None' for them which will get all.
846 You can also pass -1 for stop which will go to the end.
847 Returns the particles
850 kwds[
'resolution'] = resolution
851 if type(tuple_selection)
is str:
852 kwds[
'molecule'] = tuple_selection
853 elif type(tuple_selection)
is tuple:
854 rbegin = tuple_selection[0]
855 rend = tuple_selection[1]
856 kwds[
'molecule'] = tuple_selection[2]
858 copynum = tuple_selection[3]
859 if copynum
is not None:
860 kwds[
'copy_index'] = copynum
864 statenum = tuple_selection[4]
865 if statenum
is not None:
866 kwds[
'state_index'] = statenum
873 residue_indexes=range(1,rbegin),
875 return s.get_selected_particles()
877 kwds[
'residue_indexes'] = range(rbegin,rend+1)
879 return s.get_selected_particles()
883 def get_db_from_csv(csvfilename):
886 with open(csvfilename)
as fh:
887 csvr = csv.DictReader(fh)
894 """Store the representations for a system."""
899 self.root_hierarchy_dict = {}
900 self.preroot_fragment_hierarchy_dict = {}
901 self.particle_to_name = {}
904 def add_name(self, name):
905 if name
not in self.db:
908 def add_residue_number(self, name, resn):
911 if resn
not in self.db[name]:
912 self.db[name][resn] = {}
914 def add_resolution(self, name, resn, resolution):
916 resolution = float(resolution)
918 self.add_residue_number(name, resn)
919 if resolution
not in self.db[name][resn]:
920 self.db[name][resn][resolution] = []
924 resolution = float(resolution)
926 self.add_residue_number(name, resn)
927 self.add_resolution(name, resn, resolution)
928 self.db[name][resn][resolution] += particles
930 (rh, prf) = self.get_root_hierarchy(p)
931 self.root_hierarchy_dict[p] = rh
932 self.preroot_fragment_hierarchy_dict[p] = prf
933 self.particle_to_name[p] = name
934 if self.model
is None:
935 self.model = particles[0].get_model()
941 names = list(self.db.keys())
947 resolution = float(resolution)
948 return self.db[name][resn][resolution]
950 def get_particles_at_closest_resolution(self, name, resn, resolution):
952 resolution = float(resolution)
953 closestres = min(self.get_residue_resolutions(name, resn),
954 key=
lambda x: abs(float(x) - float(resolution)))
955 return self.get_particles(name, resn, closestres)
957 def get_residue_resolutions(self, name, resn):
959 resolutions = list(self.db[name][resn].keys())
963 def get_molecule_resolutions(self, name):
965 for resn
in self.db[name]:
966 resolutions.update(list(self.db[name][resn].keys()))
970 def get_residue_numbers(self, name):
971 residue_numbers = list(self.db[name].keys())
972 residue_numbers.sort()
973 return residue_numbers
975 def get_particles_by_resolution(self, name, resolution):
976 resolution = float(resolution)
978 for resn
in self.get_residue_numbers(name):
979 result = self.get_particles_at_closest_resolution(
983 pstemp = [p
for p
in result
if p
not in particles]
987 def get_all_particles_by_resolution(self, resolution):
988 resolution = float(resolution)
990 for name
in self.get_names():
991 particles += self.get_particles_by_resolution(name, resolution)
994 def get_root_hierarchy(self, particle):
995 prerootfragment = particle
1005 prerootfragment = particle
1011 def get_all_root_hierarchies_by_resolution(self, resolution):
1013 resolution = float(resolution)
1014 particles = self.get_all_particles_by_resolution(resolution)
1016 rh = self.root_hierarchy_dict[p]
1017 if rh
not in hierarchies:
1021 def get_preroot_fragments_by_resolution(self, name, resolution):
1023 resolution = float(resolution)
1024 particles = self.get_particles_by_resolution(name, resolution)
1026 fr = self.preroot_fragment_hierarchy_dict[p]
1027 if fr
not in fragments:
1028 fragments.append(fr)
1031 def show(self, name):
1033 for resn
in self.get_residue_numbers(name):
1035 for resolution
in self.get_residue_resolutions(name, resn):
1036 print(
"----", resolution)
1037 for p
in self.get_particles(name, resn, resolution):
1038 print(
"--------", p.get_name())
1042 '''Return the component name provided a particle and a list of names'''
1044 protname = root.get_name()
1046 while not protname
in list_of_names:
1047 root0 = root.get_parent()
1050 protname = root0.get_name()
1055 if "Beads" in protname:
1058 return (protname, is_a_bead)
1063 Retrieve the residue indexes for the given particle.
1065 The particle must be an instance of Fragment,Residue, Atom or Molecule
1066 or else returns an empty list
1077 resind_tmp=IMP.pmi1.tools.OrderedSet()
1083 resind=list(resind_tmp)
1089 def sort_by_residues(particles):
1092 sorted_particles_residues = sorted(
1094 key=
lambda tup: tup[1])
1095 particles = [p[0]
for p
in sorted_particles_residues]
1099 def get_residue_to_particle_map(particles):
1101 particles = sort_by_residues(particles)
1104 return dict(zip(particles_residues, particles))
1112 """Synchronize data over a parallel run"""
1113 from mpi4py
import MPI
1114 comm = MPI.COMM_WORLD
1115 rank = comm.Get_rank()
1116 number_of_processes = comm.size
1119 comm.send(data, dest=0, tag=11)
1122 for i
in range(1, number_of_processes):
1123 data_tmp = comm.recv(source=i, tag=11)
1124 if type(data) == list:
1126 elif type(data) == dict:
1127 data.update(data_tmp)
1129 raise TypeError(
"data not supported, use list or dictionaries")
1131 for i
in range(1, number_of_processes):
1132 comm.send(data, dest=i, tag=11)
1135 data = comm.recv(source=0, tag=11)
1139 """Synchronize data over a parallel run"""
1140 from mpi4py
import MPI
1141 comm = MPI.COMM_WORLD
1142 rank = comm.Get_rank()
1143 number_of_processes = comm.size
1146 comm.send(data, dest=0, tag=11)
1148 for i
in range(1, number_of_processes):
1149 data_tmp = comm.recv(source=i, tag=11)
1151 data[k]+=data_tmp[k]
1153 for i
in range(1, number_of_processes):
1154 comm.send(data, dest=i, tag=11)
1157 data = comm.recv(source=0, tag=11)
1167 Yield all sublists of length >= lmin and <= lmax
1175 for j
in range(i + 1, n):
1176 if len(l[i:j]) <= lmax
and len(l[i:j]) >= lmin:
1180 def flatten_list(l):
1181 return [item
for sublist
in l
for item
in sublist]
1185 """ Yield successive length-sized chunks from a list.
1187 for i
in range(0, len(list), length):
1188 yield list[i:i + length]
1191 def chunk_list_into_segments(seq, num):
1193 avg = len(seq) / float(num)
1197 while last < len(seq):
1198 out.append(seq[int(last):int(last + avg)])
1206 ''' This class stores integers
1207 in ordered compact lists eg:
1209 the methods help splitting and merging the internal lists
1211 s=Segments([1,2,3]) is [[1,2,3]]
1212 s.add(4) is [[1,2,3,4]] (add right)
1213 s.add(3) is [[1,2,3,4]] (item already existing)
1214 s.add(7) is [[1,2,3,4],[7]] (new list)
1215 s.add([8,9]) is [[1,2,3,4],[7,8,9]] (add item right)
1216 s.add([5,6]) is [[1,2,3,4,5,6,7,8,9]] (merge)
1217 s.remove(3) is [[1,2],[4,5,6,7,8,9]] (split)
1222 '''index can be a integer or a list of integers '''
1223 if type(index)
is int:
1225 elif type(index)
is list:
1226 self.segs=[[index[0]]]
1231 '''index can be a integer or a list of integers '''
1232 if type(index)
is int:
1235 for n,s
in enumerate(self.segs):
1243 if mergeright
is None and mergeleft
is None:
1244 self.segs.append([index])
1245 if not mergeright
is None and mergeleft
is None:
1246 self.segs[mergeright].append(index)
1247 if not mergeleft
is None and mergeright
is None:
1248 self.segs[mergeleft]=[index]+self.segs[mergeleft]
1249 if not mergeleft
is None and not mergeright
is None:
1250 self.segs[mergeright]=self.segs[mergeright]+[index]+self.segs[mergeleft]
1251 del self.segs[mergeleft]
1253 for n
in range(len(self.segs)):
1256 self.segs.sort(key=
lambda tup: tup[0])
1258 elif type(index)
is list:
1263 '''index can be a integer'''
1264 for n,s
in enumerate(self.segs):
1271 i=self.segs[n].index(index)
1273 self.segs.append(s[i+1:])
1274 for n
in range(len(self.segs)):
1276 if len(self.segs[n])==0:
1278 self.segs.sort(key=
lambda tup: tup[0])
1281 ''' Returns a flatten list '''
1282 return [item
for sublist
in self.segs
for item
in sublist]
1286 for seg
in self.segs:
1287 ret_tmp+=str(seg[0])+
"-"+str(seg[-1])+
","
1288 ret=ret_tmp[:-1]+
"]"
1300 Apply a translation to a hierarchy along the input vector.
1304 if type(translation_vector) == list:
1323 def translate_hierarchies(hierarchies, translation_vector):
1324 for h
in hierarchies:
1328 def translate_hierarchies_to_reference_frame(hierarchies):
1333 for h
in hierarchies:
1343 IMP.pmi1.tools.translate_hierarchies(hierarchies, (-xc, -yc, -zc))
1350 def normal_density_function(expected_value, sigma, x):
1352 1 / math.sqrt(2 * math.pi) / sigma *
1353 math.exp(-(x - expected_value) ** 2 / 2 / sigma / sigma)
1357 def log_normal_density_function(expected_value, sigma, x):
1359 1 / math.sqrt(2 * math.pi) / sigma / x *
1360 math.exp(-(math.log(x / expected_value) ** 2 / 2 / sigma / sigma))
1364 def get_random_residue_pairs(representation, resolution,
1367 avoid_same_particles=
False,
1372 names=list(representation.hier_dict.keys())
1375 prot = representation.hier_dict[name]
1376 particles +=
select(representation,name=name,resolution=resolution)
1377 random_residue_pairs = []
1378 while len(random_residue_pairs)<=number:
1379 p1 = random.choice(particles)
1380 p2 = random.choice(particles)
1381 if max_distance
is not None and \
1386 if r1==r2
and avoid_same_particles:
continue
1387 name1 = representation.get_prot_name_from_particle(p1)
1388 name2 = representation.get_prot_name_from_particle(p2)
1389 random_residue_pairs.append((name1, r1, name2, r2))
1391 return random_residue_pairs
1394 def get_random_data_point(
1400 begin_end_nbins_tuple,
1405 begin = begin_end_nbins_tuple[0]
1406 end = begin_end_nbins_tuple[1]
1407 nbins = begin_end_nbins_tuple[2]
1410 fmod_grid =
get_grid(begin, end, nbins,
True)
1412 fmod_grid = get_log_grid(begin, end, nbins)
1419 for i
in range(0, ntrials):
1420 a.append([random.random(),
True])
1423 for j
in range(1, len(fmod_grid)):
1425 fjm1 = fmod_grid[j - 1]
1429 pj = normal_density_function(expected_value, sigma, fj)
1430 pjm1 = normal_density_function(expected_value, sigma, fjm1)
1432 pj = log_normal_density_function(expected_value, sigma, fj)
1433 pjm1 = log_normal_density_function(expected_value, sigma, fjm1)
1435 norm += (pj + pjm1) / 2.0 * df
1441 for i
in range(len(cumul)):
1444 if (aa[0] <= cumul[i] / norm
and aa[1]):
1445 random_points.append(
1446 int(fmod_grid[i] / sensitivity) * sensitivity)
1450 random_points = [expected_value] * ntrials
1452 for i
in range(len(random_points)):
1453 if random.random() < outlierprob:
1454 a = random.uniform(begin, end)
1455 random_points[i] = int(a / sensitivity) * sensitivity
1456 print(random_points)
1458 for i in range(ntrials):
1459 if random.random() > OUTLIERPROB_:
1460 r=truncnorm.rvs(0.0,1.0,expected_value,BETA_)
1461 if r>1.0: print r,expected_value,BETA_
1464 random_points.append(int(r/sensitivity)*sensitivity)
1469 for r
in random_points:
1473 rmean /= float(ntrials)
1474 rmean2 /= float(ntrials)
1475 stddev = math.sqrt(max(rmean2 - rmean * rmean, 0.))
1476 return rmean, stddev
1478 def print_multicolumn(list_of_strings, ncolumns=2, truncate=40):
1484 for i
in range(len(l) % cols):
1487 split = [l[i:i + len(l) / cols]
for i
in range(0, len(l), len(l) / cols)]
1488 for row
in zip(*split):
1489 print(
"".join(str.ljust(i, truncate)
for i
in row))
1492 '''Change color code to hexadecimal to rgb'''
1494 self._NUMERALS =
'0123456789abcdefABCDEF'
1495 self._HEXDEC = dict((v, int(v, 16))
for v
in (x+y
for x
in self._NUMERALS
for y
in self._NUMERALS))
1496 self.LOWERCASE, self.UPPERCASE =
'x',
'X'
1498 def rgb(self,triplet):
1499 return float(self._HEXDEC[triplet[0:2]]), float(self._HEXDEC[triplet[2:4]]), float(self._HEXDEC[triplet[4:6]])
1501 def triplet(self,rgb, lettercase=None):
1502 if lettercase
is None: lettercase=self.LOWERCASE
1503 return format(rgb[0]<<16 | rgb[1]<<8 | rgb[2],
'06'+lettercase)
1507 class OrderedSet(collections.MutableSet):
1509 def __init__(self, iterable=None):
1511 end += [
None, end, end]
1513 if iterable
is not None:
1517 return len(self.map)
1519 def __contains__(self, key):
1520 return key
in self.map
1523 if key
not in self.map:
1526 curr[2] = end[1] = self.map[key] = [key, curr, end]
1528 def discard(self, key):
1530 key, prev, next = self.map.pop(key)
1537 while curr
is not end:
1541 def __reversed__(self):
1544 while curr
is not end:
1548 def pop(self, last=True):
1550 raise KeyError(
'set is empty')
1552 key = self.end[1][0]
1554 key = self.end[2][0]
1560 return '%s()' % (self.__class__.__name__,)
1561 return '%s(%r)' % (self.__class__.__name__, list(self))
1563 def __eq__(self, other):
1564 if isinstance(other, OrderedSet):
1565 return len(self) == len(other)
and list(self) == list(other)
1566 return set(self) == set(other)
1570 """Store objects in order they were added, but with default type.
1571 Source: http://stackoverflow.com/a/4127426/2608793
1573 def __init__(self, *args, **kwargs):
1575 self.default_factory =
None
1577 if not (args[0]
is None or callable(args[0])):
1578 raise TypeError(
'first argument must be callable or None')
1579 self.default_factory = args[0]
1581 super(OrderedDefaultDict, self).__init__(*args, **kwargs)
1583 def __missing__ (self, key):
1584 if self.default_factory
is None:
1586 self[key] = default = self.default_factory()
1589 def __reduce__(self):
1590 args = (self.default_factory,)
if self.default_factory
else ()
1591 return self.__class__, args,
None,
None, self.iteritems()
1596 """Extract frame from RMF file and fill coordinates. Must be identical topology.
1597 @param hier The (System) hierarchy to fill (e.g. after you've built it)
1598 @param rmf_fn The file to extract from
1599 @param frame_num The frame number to extract
1601 rh = RMF.open_rmf_file_read_only(rmf_fn)
1609 selection_tuple=
None,
1610 warn_about_slices=
True):
1611 """Adapt things for PMI (degrees of freedom, restraints, ...)
1612 Returns list of list of hierarchies, separated into Molecules if possible.
1613 The input can be a list, or a list of lists (iterable of ^1 or iterable of ^2)
1614 (iterable of ^2) Hierarchy -> returns input as list of list of hierarchies,
1615 only one entry, not grouped by molecules.
1616 (iterable of ^2) PMI::System/State/Molecule/TempResidue ->
1617 returns residue hierarchies, grouped in molecules, at requested resolution
1618 @param stuff Can be one of the following inputs:
1619 IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
1620 Must be uniform input, however. No mixing object types.
1621 @param pmi_resolution For selecting, only does it if you pass PMI objects. Set it to "all"
1622 if you want all resolutions!
1623 @param flatten Set to True if you just want all hierarchies in one list.
1624 @param warn_about_slices Print a warning if you are requesting only part of a bead.
1625 Sometimes you just don't care!
1626 \note since this relies on IMP::atom::Selection, this will not return any objects if they weren't built!
1627 But there should be no problem if you request unbuilt residues, they should be ignored.
1633 if hasattr(stuff,
'__iter__'):
1639 if all(hasattr(el,
'__iter__')
for el
in thelist):
1640 thelist = [i
for sublist
in thelist
for i
in sublist]
1641 elif any(hasattr(el,
'__iter__')
for el
in thelist):
1642 raise Exception(
'input_adaptor: input_object must be a list or a list of lists')
1651 except NotImplementedError:
1657 is_temp_residue=
False
1662 if is_system
or is_state
or is_molecule
or is_temp_residue:
1667 for system
in stuff:
1668 for state
in system.get_states():
1669 mdict = state.get_molecules()
1670 for molname
in mdict:
1671 for copy
in mdict[molname]:
1672 indexes_per_mol[copy] += [r.get_index()
for r
in copy.get_residues()]
1675 mdict = state.get_molecules()
1676 for molname
in mdict:
1677 for copy
in mdict[molname]:
1678 indexes_per_mol[copy] += [r.get_index()
for r
in copy.get_residues()]
1680 for molecule
in stuff:
1681 indexes_per_mol[molecule] += [r.get_index()
for r
in molecule.get_residues()]
1682 elif is_temp_residue:
1683 for tempres
in stuff:
1684 indexes_per_mol[tempres.get_molecule()].append(tempres.get_index())
1685 for mol
in indexes_per_mol:
1686 if pmi_resolution==
'all':
1690 residue_indexes=indexes_per_mol[mol])
1693 resolution=pmi_resolution,
1694 residue_indexes=indexes_per_mol[mol])
1695 ps = sel.get_selected_particles()
1698 if warn_about_slices:
1699 rset = set(indexes_per_mol[mol])
1703 if not fset <= rset:
1710 resbreak = maxf
if minf==minset
else minset-1
1711 print(
'WARNING: You are trying to select only part of the bead %s:%i-%i.\n'
1712 'The residues you requested are %i-%i. You can fix this by:\n'
1713 '1) requesting the whole bead/none of it or\n'
1714 '2) break the bead up by passing bead_extra_breaks=[\'%i\'] in '
1715 'molecule.add_representation()'
1716 %(mol.get_name(),minset,maxset,minf,maxf,resbreak))
1721 if pmi_resolution==
'all':
1729 hier_list = [hier_list]
1731 raise Exception(
'input_adaptor: you passed something of wrong type or a list with mixed types')
1733 if flatten
and pmi_input:
1734 return [h
for sublist
in hier_list
for h
in sublist]
1740 """Returns sequence-sorted segments array, each containing the first particle
1741 the last particle and the first residue index."""
1743 from operator
import itemgetter
1746 raise Exception(
"IMP.pmi1.tools.get_sorted_segments: only pass stuff from one Molecule, please")
1762 SortedSegments.append((start, end, startres))
1763 SortedSegments = sorted(SortedSegments, key=itemgetter(2))
1764 return SortedSegments
1767 """Decorate the sequence-consecutive particles from a PMI2 molecule with a bond,
1768 so that they appear connected in the rmf file"""
1770 for x
in range(len(SortedSegments) - 1):
1772 last = SortedSegments[x][1]
1773 first = SortedSegments[x + 1][0]
1775 p1 = last.get_particle()
1776 p2 = first.get_particle()
1789 """This class converts three to one letter codes, and return X for any unknown codes"""
1790 def __init__(self,is_nucleic=False):
1793 threetoone = {
'ALA':
'A',
'ARG':
'R', 'ASN': 'N', 'ASP': 'D',
1794 'CYS':
'C',
'GLU':
'E',
'GLN':
'Q',
'GLY':
'G',
1795 'HIS':
'H',
'ILE':
'I',
'LEU':
'L',
'LYS':
'K',
1796 'MET':
'M',
'PHE':
'F',
'PRO':
'P',
'SER':
'S',
1797 'THR':
'T',
'TRP':
'W',
'TYR':
'Y',
'VAL':
'V',
'UNK':
'X'}
1799 threetoone = {
'ADE':
'A',
'URA':
'U', 'CYT': 'C', 'GUA': 'G',
1800 'THY':
'T',
'UNK':
'X'}
1802 defaultdict.__init__(self,
lambda:
"X", threetoone)
1807 def get_residue_type_from_one_letter_code(code,is_nucleic=False):
1810 for k
in threetoone:
1811 one_to_three[threetoone[k]] = k
1816 """ Just get the leaves from a list of hierarchies """
1817 lvs = list(itertools.chain.from_iterable(
IMP.atom.get_leaves(item)
for item
in list_of_hs))
1824 """Perform selection using the usual keywords but return ALL resolutions (BEADS and GAUSSIANS).
1825 Returns in flat list!
1830 if hier
is not None:
1833 print(
"WARNING: You passed nothing to select_at_all_resolutions()")
1840 raise Exception(
'select_at_all_resolutions: you have to pass an IMP Hierarchy')
1842 raise Exception(
'select_at_all_resolutions: you have to pass an IMP Hierarchy')
1843 if 'resolution' in kwargs
or 'representation_type' in kwargs:
1844 raise Exception(
"don't pass resolution or representation_type to this function")
1846 representation_type=IMP.atom.BALLS,**kwargs)
1848 representation_type=IMP.atom.DENSITIES,**kwargs)
1849 ret |= OrderedSet(selB.get_selected_particles())
1850 ret |= OrderedSet(selD.get_selected_particles())
1859 """Utility to retrieve particles from a hierarchy within a
1860 zone around a set of ps.
1861 @param hier The hierarchy in which to look for neighbors
1862 @param target_ps The particles for zoning
1863 @param sel_zone The maximum distance
1864 @param entire_residues If True, will grab entire residues
1865 @param exclude_backbone If True, will only return sidechain particles
1869 backbone_types=[
'C',
'N',
'CB',
'O']
1870 if exclude_backbone:
1872 for n
in backbone_types])
1873 test_ps = test_sel.get_selected_particles()
1874 nn = IMP.algebra.NearestNeighbor3D([
IMP.core.XYZ(p).get_coordinates()
1877 for target
in target_ps:
1878 zone|=set(nn.get_in_ball(
IMP.core.XYZ(target).get_coordinates(),sel_zone))
1879 zone_ps = [test_ps[z]
for z
in zone]
1884 zone_ps = [h.get_particle()
for h
in final_ps]
1889 """Returns unique objects in original order"""
1893 if not hasattr(hiers,
'__iter__'):
1900 rbs_ordered.append(rb)
1905 rbs_ordered.append(rb)
1909 return rbs_ordered,beads
1912 "This function returns the parent molecule hierarchies of given objects"
1913 stuff=
input_adaptor(input_objects, pmi_resolution=
'all',flatten=
True)
1918 while not (is_molecule
or is_root):
1919 root=IMP.atom.get_root(h)
1926 return list(molecules)
1928 def get_molecules_dictionary(input_objects):
1929 moldict=defaultdict(list)
1932 moldict[name].append(mol)
1938 def get_molecules_dictionary_by_copy(input_objects):
1939 moldict=defaultdict(dict)
1943 moldict[name][c]=mol
1946 def get_selections_dictionary(input_objects):
1947 moldict=IMP.pmi1.tools.get_molecules_dictionary(input_objects)
1948 seldict=defaultdict(list)
1949 for name, mols
in moldict.items():
1955 """Given a list of PMI objects, returns all density hierarchies within
1956 these objects. The output of this function can be inputted into
1957 things such as EM restraints. This function is intended to gather density particles
1958 appended to molecules (and not other hierarchies which might have been appended to the root node directly).
1965 densities+=
IMP.atom.Selection(i,representation_type=IMP.atom.DENSITIES).get_selected_particles()
1969 max_translation=300., max_rotation=2.0 * pi,
1970 avoidcollision_rb=
True, avoidcollision_fb=
False,
1971 cutoff=10.0, niterations=100,
1973 excluded_rigid_bodies=[],
1974 hierarchies_excluded_from_collision=[],
1975 hierarchies_included_in_collision=[],
1977 return_debug=
False):
1978 """Shuffle particles. Used to restart the optimization.
1979 The configuration of the system is initialized by placing each
1980 rigid body and each bead randomly in a box. If `bounding_box` is
1981 specified, the particles are placed inside this box; otherwise, each
1982 particle is displaced by up to max_translation angstroms, and randomly
1983 rotated. Effort is made to place particles far enough from each other to
1984 prevent any steric clashes.
1985 @param objects Can be one of the following inputs:
1986 IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set of them
1987 @param max_translation Max translation (rbs and flexible beads)
1988 @param max_rotation Max rotation (rbs only)
1989 @param avoidcollision_rb check if the particle/rigid body was
1990 placed close to another particle; uses the optional
1991 arguments cutoff and niterations
1992 @param avoidcollision_fb Advanced. Generally you want this False because it's hard to shuffle beads.
1993 @param cutoff Distance less than this is a collision
1994 @param niterations How many times to try avoiding collision
1995 @param bounding_box Only shuffle particles within this box. Defined by ((x1,y1,z1),(x2,y2,z2)).
1996 @param excluded_rigid_bodies Don't shuffle these rigid body objects
1997 @param hierarchies_excluded_from_collision Don't count collision with these bodies
1998 @param hierarchies_included_in_collision Hierarchies that are not shuffled, but should be included in collision calculation (for fixed regions)
1999 @param verbose Give more output
2000 \note Best to only call this function after you've set up degrees of freedom
2001 For debugging purposes, returns: <shuffled indexes>, <collision avoided indexes>
2006 pmi_resolution=
'all',
2009 if len(rigid_bodies)>0:
2010 mdl = rigid_bodies[0].get_model()
2011 elif len(flexible_beads)>0:
2012 mdl = flexible_beads[0].get_model()
2014 raise Exception(
"Could not find any particles in the hierarchy")
2015 if len(rigid_bodies) == 0:
2016 print(
"shuffle_configuration: rigid bodies were not intialized")
2020 gcpf.set_distance(cutoff)
2024 pmi_resolution=
'all',
2028 pmi_resolution=
'all',
2031 collision_excluded_idxs = set([l.get_particle().
get_index()
for h
in collision_excluded_hierarchies \
2034 collision_included_idxs = set([l.get_particle().
get_index()
for h
in collision_included_hierarchies \
2041 all_idxs.append(p.get_particle_index())
2043 collision_excluded_idxs.add(p.get_particle_index())
2045 if bounding_box
is not None:
2046 ((x1, y1, z1), (x2, y2, z2)) = bounding_box
2051 all_idxs = set(all_idxs) | collision_included_idxs
2052 all_idxs = all_idxs - collision_excluded_idxs
2054 print(
'shuffling', len(rigid_bodies),
'rigid bodies')
2055 for rb
in rigid_bodies:
2056 if rb
not in excluded_rigid_bodies:
2058 if avoidcollision_rb:
2059 rb_idxs = set(rb.get_member_particle_indexes()) - \
2060 collision_excluded_idxs
2061 other_idxs = all_idxs - rb_idxs
2067 while niter < niterations:
2068 rbxyz = (rb.get_x(), rb.get_y(), rb.get_z())
2088 debug.append([rb, other_idxs
if avoidcollision_rb
else set()])
2092 if avoidcollision_rb:
2094 npairs = len(gcpf.get_close_pairs(mdl,
2103 print(
"shuffle_configuration: rigid body placed close to other %d particles, trying again..." % npairs)
2104 print(
"shuffle_configuration: rigid body name: " + rb.get_name())
2105 if niter == niterations:
2106 raise ValueError(
"tried the maximum number of iterations to avoid collisions, increase the distance cutoff")
2110 print(
'shuffling', len(flexible_beads),
'flexible beads')
2111 for fb
in flexible_beads:
2113 if avoidcollision_fb:
2115 other_idxs = all_idxs - fb_idxs
2121 while niter < niterations:
2134 xyz = memb.get_internal_coordinates()
2139 rf = memb.get_rigid_body().get_reference_frame()
2140 glob_to_int = rf.get_transformation_from()
2141 memb.set_internal_coordinates(
2142 glob_to_int.get_transformed(translation))
2144 xyz_transformed=transformation.get_transformed(xyz)
2145 memb.set_internal_coordinates(xyz_transformed)
2146 debug.append([xyz,other_idxs
if avoidcollision_fb
else set()])
2153 debug.append([d,other_idxs
if avoidcollision_fb
else set()])
2159 if avoidcollision_fb:
2161 npairs = len(gcpf.get_close_pairs(mdl,
2169 print(
"shuffle_configuration: floppy body placed close to other %d particles, trying again..." % npairs)
2170 if niter == niterations:
2171 raise ValueError(
"tried the maximum number of iterations to avoid collisions, increase the distance cutoff")
2177 class ColorHierarchy(object):
2179 def __init__(self,hier):
2180 import matplotlib
as mpl
2182 import matplotlib.pyplot
as plt
2186 hier.ColorHierarchy=self
2190 self.method=self.nochange
2198 def get_color(self,fl):
2201 def get_log_scale(self,fl):
2204 return math.log(fl+eps)
2206 def color_by_resid(self):
2207 self.method=self.color_by_resid
2208 self.scheme=self.mpl.cm.rainbow
2209 for mol
in self.mols:
2215 c=self.get_color(float(ri)/self.last)
2219 avr=sum(ris)/len(ris)
2220 c=self.get_color(float(avr)/self.last)
2223 def color_by_uncertainty(self):
2224 self.method=self.color_by_uncertainty
2225 self.scheme=self.mpl.cm.jet
2232 self.first=self.get_log_scale(1.0)
2233 self.last=self.get_log_scale(100.0)
2235 value=self.get_log_scale(unc_dict[p])
2236 if value>=self.last: value=self.last
2237 if value<=self.first: value=self.first
2238 c=self.get_color((value-self.first)/(self.last-self.first))
2241 def get_color_bar(self,filename):
2242 import matplotlib
as mpl
2244 import matplotlib.pyplot
as plt
2247 fig = plt.figure(figsize=(8, 3))
2248 ax1 = fig.add_axes([0.05, 0.80, 0.9, 0.15])
2251 norm = mpl.colors.Normalize(vmin=0.0, vmax=1.0)
2253 if self.method == self.color_by_uncertainty:
2254 angticks=[1.0,2.5,5.0,10.0,25.0,50.0,100.0]
2258 vvalue=(self.get_log_scale(at)-self.first)/(self.last-self.first)
2259 if vvalue <= 1.0
and vvalue >= 0.0:
2260 vvalues.append(vvalue)
2261 marks.append(str(at))
2262 cb1 = mpl.colorbar.ColorbarBase(ax1, cmap=cmap,
2265 orientation=
'horizontal')
2266 print(self.first,self.last,marks,vvalues)
2267 cb1.ax.set_xticklabels(marks)
2268 cb1.set_label(
'Angstorm')
2269 plt.savefig(filename, dpi=150, transparent=
True)
2276 """Given a chimera color name, return RGB"""
2277 d = {
'aquamarine': (0.4980392156862745, 1.0, 0.8313725490196079),
2278 'black': (0.0, 0.0, 0.0),
2279 'blue': (0.0, 0.0, 1.0),
2280 'brown': (0.6470588235294118, 0.16470588235294117, 0.16470588235294117),
2281 'chartreuse': (0.4980392156862745, 1.0, 0.0),
2282 'coral': (1.0, 0.4980392156862745, 0.3137254901960784),
2283 'cornflower blue': (0.39215686274509803, 0.5843137254901961, 0.9294117647058824),
2284 'cyan': (0.0, 1.0, 1.0),
2285 'dark cyan': (0.0, 0.5450980392156862, 0.5450980392156862),
2286 'dark gray': (0.6627450980392157, 0.6627450980392157, 0.6627450980392157),
2287 'dark green': (0.0, 0.39215686274509803, 0.0),
2288 'dark khaki': (0.7411764705882353, 0.7176470588235294, 0.4196078431372549),
2289 'dark magenta': (0.5450980392156862, 0.0, 0.5450980392156862),
2290 'dark olive green': (0.3333333333333333, 0.4196078431372549, 0.1843137254901961),
2291 'dark red': (0.5450980392156862, 0.0, 0.0),
2292 'dark slate blue': (0.2823529411764706, 0.23921568627450981, 0.5450980392156862),
2293 'dark slate gray': (0.1843137254901961, 0.30980392156862746, 0.30980392156862746),
2294 'deep pink': (1.0, 0.0784313725490196, 0.5764705882352941),
2295 'deep sky blue': (0.0, 0.7490196078431373, 1.0),
2296 'dim gray': (0.4117647058823529, 0.4117647058823529, 0.4117647058823529),
2297 'dodger blue': (0.11764705882352941, 0.5647058823529412, 1.0),
2298 'firebrick': (0.6980392156862745, 0.13333333333333333, 0.13333333333333333),
2299 'forest green': (0.13333333333333333, 0.5450980392156862, 0.13333333333333333),
2300 'gold': (1.0, 0.8431372549019608, 0.0),
2301 'goldenrod': (0.8549019607843137, 0.6470588235294118, 0.12549019607843137),
2302 'gray': (0.7450980392156863, 0.7450980392156863, 0.7450980392156863),
2303 'green': (0.0, 1.0, 0.0),
2304 'hot pink': (1.0, 0.4117647058823529, 0.7058823529411765),
2305 'khaki': (0.9411764705882353, 0.9019607843137255, 0.5490196078431373),
2306 'light blue': (0.6784313725490196, 0.8470588235294118, 0.9019607843137255),
2307 'light gray': (0.8274509803921568, 0.8274509803921568, 0.8274509803921568),
2308 'light green': (0.5647058823529412, 0.9333333333333333, 0.5647058823529412),
2309 'light sea green': (0.12549019607843137, 0.6980392156862745, 0.6666666666666666),
2310 'lime green': (0.19607843137254902, 0.803921568627451, 0.19607843137254902),
2311 'magenta': (1.0, 0.0, 1.0),
2312 'medium blue': (0.19607843137254902, 0.19607843137254902, 0.803921568627451),
2313 'medium purple': (0.5764705882352941, 0.4392156862745098, 0.8588235294117647),
2314 'navy blue': (0.0, 0.0, 0.5019607843137255),
2315 'olive drab': (0.4196078431372549, 0.5568627450980392, 0.13725490196078433),
2316 'orange red': (1.0, 0.27058823529411763, 0.0),
2317 'orange': (1.0, 0.4980392156862745, 0.0),
2318 'orchid': (0.8549019607843137, 0.4392156862745098, 0.8392156862745098),
2319 'pink': (1.0, 0.7529411764705882, 0.796078431372549),
2320 'plum': (0.8666666666666667, 0.6274509803921569, 0.8666666666666667),
2321 'purple': (0.6274509803921569, 0.12549019607843137, 0.9411764705882353),
2322 'red': (1.0, 0.0, 0.0),
2323 'rosy brown': (0.7372549019607844, 0.5607843137254902, 0.5607843137254902),
2324 'salmon': (0.9803921568627451, 0.5019607843137255, 0.4470588235294118),
2325 'sandy brown': (0.9568627450980393, 0.6431372549019608, 0.3764705882352941),
2326 'sea green': (0.1803921568627451, 0.5450980392156862, 0.3411764705882353),
2327 'sienna': (0.6274509803921569, 0.3215686274509804, 0.17647058823529413),
2328 'sky blue': (0.5294117647058824, 0.807843137254902, 0.9215686274509803),
2329 'slate gray': (0.4392156862745098, 0.5019607843137255, 0.5647058823529412),
2330 'spring green': (0.0, 1.0, 0.4980392156862745),
2331 'steel blue': (0.27450980392156865, 0.5098039215686274, 0.7058823529411765),
2332 'tan': (0.8235294117647058, 0.7058823529411765, 0.5490196078431373),
2333 'turquoise': (0.25098039215686274, 0.8784313725490196, 0.8156862745098039),
2334 'violet red': (0.8156862745098039, 0.12549019607843137, 0.5647058823529412),
2335 'white': (1.0, 1.0, 1.0),
2336 'yellow': (1.0, 1.0, 0.0)}
2342 "reds":[(
"maroon",
"#800000",(128,0,0)),(
"dark red",
"#8B0000",(139,0,0)),
2343 (
"brown",
"#A52A2A",(165,42,42)),(
"firebrick",
"#B22222",(178,34,34)),
2344 (
"crimson",
"#DC143C",(220,20,60)),(
"red",
"#FF0000",(255,0,0)),
2345 (
"tomato",
"#FF6347",(255,99,71)),(
"coral",
"#FF7F50",(255,127,80)),
2346 (
"indian red",
"#CD5C5C",(205,92,92)),(
"light coral",
"#F08080",(240,128,128)),
2347 (
"dark salmon",
"#E9967A",(233,150,122)),(
"salmon",
"#FA8072",(250,128,114)),
2348 (
"light salmon",
"#FFA07A",(255,160,122)),(
"orange red",
"#FF4500",(255,69,0)),
2349 (
"dark orange",
"#FF8C00",(255,140,0))],
2350 "yellows":[(
"orange",
"#FFA500",(255,165,0)),(
"gold",
"#FFD700",(255,215,0)),
2351 (
"dark golden rod",
"#B8860B",(184,134,11)),(
"golden rod",
"#DAA520",(218,165,32)),
2352 (
"pale golden rod",
"#EEE8AA",(238,232,170)),(
"dark khaki",
"#BDB76B",(189,183,107)),
2353 (
"khaki",
"#F0E68C",(240,230,140)),(
"olive",
"#808000",(128,128,0)),
2354 (
"yellow",
"#FFFF00",(255,255,0)),(
"antique white",
"#FAEBD7",(250,235,215)),
2355 (
"beige",
"#F5F5DC",(245,245,220)),(
"bisque",
"#FFE4C4",(255,228,196)),
2356 (
"blanched almond",
"#FFEBCD",(255,235,205)),(
"wheat",
"#F5DEB3",(245,222,179)),
2357 (
"corn silk",
"#FFF8DC",(255,248,220)),(
"lemon chiffon",
"#FFFACD",(255,250,205)),
2358 (
"light golden rod yellow",
"#FAFAD2",(250,250,210)),(
"light yellow",
"#FFFFE0",(255,255,224))],
2359 "greens":[(
"yellow green",
"#9ACD32",(154,205,50)),(
"dark olive green",
"#556B2F",(85,107,47)),
2360 (
"olive drab",
"#6B8E23",(107,142,35)),(
"lawn green",
"#7CFC00",(124,252,0)),
2361 (
"chart reuse",
"#7FFF00",(127,255,0)),(
"green yellow",
"#ADFF2F",(173,255,47)),
2362 (
"dark green",
"#006400",(0,100,0)),(
"green",
"#008000",(0,128,0)),
2363 (
"forest green",
"#228B22",(34,139,34)),(
"lime",
"#00FF00",(0,255,0)),
2364 (
"lime green",
"#32CD32",(50,205,50)),(
"light green",
"#90EE90",(144,238,144)),
2365 (
"pale green",
"#98FB98",(152,251,152)),(
"dark sea green",
"#8FBC8F",(143,188,143)),
2366 (
"medium spring green",
"#00FA9A",(0,250,154)),(
"spring green",
"#00FF7F",(0,255,127)),
2367 (
"sea green",
"#2E8B57",(46,139,87)),(
"medium aqua marine",
"#66CDAA",(102,205,170)),
2368 (
"medium sea green",
"#3CB371",(60,179,113)),(
"light sea green",
"#20B2AA",(32,178,170)),
2369 (
"dark slate gray",
"#2F4F4F",(47,79,79)),(
"teal",
"#008080",(0,128,128)),
2370 (
"dark cyan",
"#008B8B",(0,139,139))],
2371 "blues":[(
"dark turquoise",
"#00CED1",(0,206,209)),
2372 (
"turquoise",
"#40E0D0",(64,224,208)),(
"medium turquoise",
"#48D1CC",(72,209,204)),
2373 (
"pale turquoise",
"#AFEEEE",(175,238,238)),(
"aqua marine",
"#7FFFD4",(127,255,212)),
2374 (
"powder blue",
"#B0E0E6",(176,224,230)),(
"cadet blue",
"#5F9EA0",(95,158,160)),
2375 (
"steel blue",
"#4682B4",(70,130,180)),(
"corn flower blue",
"#6495ED",(100,149,237)),
2376 (
"deep sky blue",
"#00BFFF",(0,191,255)),(
"dodger blue",
"#1E90FF",(30,144,255)),
2377 (
"light blue",
"#ADD8E6",(173,216,230)),(
"sky blue",
"#87CEEB",(135,206,235)),
2378 (
"light sky blue",
"#87CEFA",(135,206,250)),(
"midnight blue",
"#191970",(25,25,112)),
2379 (
"navy",
"#000080",(0,0,128)),(
"dark blue",
"#00008B",(0,0,139)),
2380 (
"medium blue",
"#0000CD",(0,0,205)),(
"blue",
"#0000FF",(0,0,255)),(
"royal blue",
"#4169E1",(65,105,225)),
2381 (
"aqua",
"#00FFFF",(0,255,255)),(
"cyan",
"#00FFFF",(0,255,255)),(
"light cyan",
"#E0FFFF",(224,255,255))],
2382 "violets":[(
"blue violet",
"#8A2BE2",(138,43,226)),(
"indigo",
"#4B0082",(75,0,130)),
2383 (
"dark slate blue",
"#483D8B",(72,61,139)),(
"slate blue",
"#6A5ACD",(106,90,205)),
2384 (
"medium slate blue",
"#7B68EE",(123,104,238)),(
"medium purple",
"#9370DB",(147,112,219)),
2385 (
"dark magenta",
"#8B008B",(139,0,139)),(
"dark violet",
"#9400D3",(148,0,211)),
2386 (
"dark orchid",
"#9932CC",(153,50,204)),(
"medium orchid",
"#BA55D3",(186,85,211)),
2387 (
"purple",
"#800080",(128,0,128)),(
"thistle",
"#D8BFD8",(216,191,216)),
2388 (
"plum",
"#DDA0DD",(221,160,221)),(
"violet",
"#EE82EE",(238,130,238)),
2389 (
"magenta / fuchsia",
"#FF00FF",(255,0,255)),(
"orchid",
"#DA70D6",(218,112,214)),
2390 (
"medium violet red",
"#C71585",(199,21,133)),(
"pale violet red",
"#DB7093",(219,112,147)),
2391 (
"deep pink",
"#FF1493",(255,20,147)),(
"hot pink",
"#FF69B4",(255,105,180)),
2392 (
"light pink",
"#FFB6C1",(255,182,193)),(
"pink",
"#FFC0CB",(255,192,203))],
2393 "browns":[(
"saddle brown",
"#8B4513",(139,69,19)),(
"sienna",
"#A0522D",(160,82,45)),
2394 (
"chocolate",
"#D2691E",(210,105,30)),(
"peru",
"#CD853F",(205,133,63)),
2395 (
"sandy brown",
"#F4A460",(244,164,96)),(
"burly wood",
"#DEB887",(222,184,135)),
2396 (
"tan",
"#D2B48C",(210,180,140)),(
"rosy brown",
"#BC8F8F",(188,143,143)),
2397 (
"moccasin",
"#FFE4B5",(255,228,181)),(
"navajo white",
"#FFDEAD",(255,222,173)),
2398 (
"peach puff",
"#FFDAB9",(255,218,185)),(
"misty rose",
"#FFE4E1",(255,228,225)),
2399 (
"lavender blush",
"#FFF0F5",(255,240,245)),(
"linen",
"#FAF0E6",(250,240,230)),
2400 (
"old lace",
"#FDF5E6",(253,245,230)),(
"papaya whip",
"#FFEFD5",(255,239,213)),
2401 (
"sea shell",
"#FFF5EE",(255,245,238))],
2402 "greys":[(
"black",
"#000000",(0,0,0)),(
"dim gray / dim grey",
"#696969",(105,105,105)),
2403 (
"gray / grey",
"#808080",(128,128,128)),(
"dark gray / dark grey",
"#A9A9A9",(169,169,169)),
2404 (
"silver",
"#C0C0C0",(192,192,192)),(
"light gray / light grey",
"#D3D3D3",(211,211,211)),
2405 (
"gainsboro",
"#DCDCDC",(220,220,220)),(
"white smoke",
"#F5F5F5",(245,245,245)),
2406 (
"white",
"#FFFFFF",(255,255,255))]}
2408 def assign_color_group(self,color_group,representation,component_names):
2409 for n,p
in enumerate(component_names):
2411 psel=s.get_selected_particles()
2412 ctuple=self.colors[color_group][n]
2413 print(
"Assigning "+p+
" to color "+ctuple[0])
2422 def get_list_distant_colors(self):
2423 cnames = [
'#F0F8FF',
'#FAEBD7',
'#00FFFF',
'#7FFFD4',
'#F0FFFF',
'#F5F5DC',
2424 '#FFE4C4',
'#000000',
'#FFEBCD',
'#0000FF',
'#8A2BE2',
'#A52A2A',
'#DEB887',
2425 '#5F9EA0',
'#7FFF00',
'#D2691E',
'#FF7F50',
'#6495ED',
'#FFF8DC',
'#DC143C',
2426 '#00FFFF',
'#00008B',
'#008B8B',
'#B8860B',
'#A9A9A9',
'#006400',
'#BDB76B',
2427 '#8B008B',
'#556B2F',
'#FF8C00',
'#9932CC',
'#8B0000',
'#E9967A',
'#8FBC8F',
2428 '#483D8B',
'#2F4F4F',
'#00CED1',
'#9400D3',
'#FF1493',
'#00BFFF',
'#696969',
2429 '#1E90FF',
'#B22222',
'#FFFAF0',
'#228B22',
'#FF00FF',
'#DCDCDC',
'#F8F8FF',
2430 '#FFD700',
'#DAA520',
'#808080',
'#008000',
'#ADFF2F',
'#F0FFF0',
'#FF69B4',
2431 '#CD5C5C',
'#4B0082',
'#FFFFF0',
'#F0E68C',
'#E6E6FA',
'#FFF0F5',
'#7CFC00',
2432 '#FFFACD',
'#ADD8E6',
'#F08080',
'#E0FFFF',
'#FAFAD2',
'#90EE90',
'#D3D3D3',
2433 '#FFB6C1',
'#FFA07A',
'#20B2AA',
'#87CEFA',
'#778899',
'#B0C4DE',
'#FFFFE0',
2434 '#00FF00',
'#32CD32',
'#FAF0E6',
'#FF00FF',
'#800000',
'#66CDAA',
'#0000CD',
2435 '#BA55D3',
'#9370DB',
'#3CB371',
'#7B68EE',
'#00FA9A',
'#48D1CC',
'#C71585',
2436 '#191970',
'#F5FFFA',
'#FFE4E1',
'#FFE4B5',
'#FFDEAD',
'#000080',
'#FDF5E6',
2437 '#808000',
'#6B8E23',
'#FFA500',
'#FF4500',
'#DA70D6',
'#EEE8AA',
'#98FB98',
2438 '#AFEEEE',
'#DB7093',
'#FFEFD5',
'#FFDAB9',
'#CD853F',
'#FFC0CB',
'#DDA0DD',
2439 '#B0E0E6',
'#800080',
'#FF0000',
'#BC8F8F',
'#4169E1',
'#8B4513',
'#FA8072',
2440 '#FAA460',
'#2E8B57',
'#FFF5EE',
'#A0522D',
'#C0C0C0',
'#87CEEB',
'#6A5ACD',
2441 '#708090',
'#FFFAFA',
'#00FF7F',
'#4682B4',
'#D2B48C',
'#008080',
'#D8BFD8',
2442 '#FF6347',
'#40E0D0',
'#EE82EE',
'#F5DEB3',
'#FFFFFF',
'#F5F5F5',
'#FFFF00',
static bool get_is_setup(const IMP::ParticleAdaptor &p)
static bool get_is_setup(const IMP::ParticleAdaptor &p)
A decorator to associate a particle with a part of a protein/DNA/RNA.
def add_script_provenance
Tag the given particle with the current Python script.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
void add_particles(RMF::FileHandle fh, const ParticlesTemp &hs)
Add resolution to a particle.
static Weight setup_particle(Model *m, ParticleIndex pi)
Set up an empty Weight.
A decorator for a particle which has bonds.
double get_drmsd(const Vector3DsOrXYZs0 &m0, const Vector3DsOrXYZs1 &m1)
Calculate distance the root mean square deviation between two sets of 3D points.
Rotation3D get_random_rotation_3d(const Rotation3D ¢er, double distance)
Pick a rotation at random near the provided one.
Legacy PMI1 module to represent, score, sample and analyze models.
Extends the functionality of IMP.atom.Molecule.
IMP::Vector< Color > Colors
Set of python classes to create a multi-state, multi-resolution IMP hierarchy.
ParticlesTemp get_particles(Model *m, const ParticleIndexes &ps)
static Surface setup_particle(Model *m, ParticleIndex pi)
void add_particle(RMF::FileHandle fh, Particle *hs)
static bool get_is_setup(const IMP::ParticleAdaptor &p)
static bool get_is_setup(const IMP::ParticleAdaptor &p)
GenericHierarchies get_leaves(Hierarchy mhd)
Get all the leaves of the bit of hierarchy.
double get_distance(XYZR a, XYZR b)
Compute the sphere distance between a and b.
def add_imp_provenance
Tag the given particle as being created by the current version of IMP.
Vector3D get_random_vector_in(const Cylinder3D &c)
Generate a random vector in a cylinder with uniform density.
Bond create_bond(Bonded a, Bonded b, Bond o)
Connect the two wrapped particles by a custom bond.
Object used to hold a set of restraints.
Add uncertainty to a particle.
static bool get_is_setup(Model *m, ParticleIndex pi)
def add_software_provenance
Tag the given particle with the software used to create it.
A decorator for keeping track of copies of a molecule.
static bool get_is_setup(Model *m, ParticleIndex pi)
ParticleIndexPairs get_indexes(const ParticlePairsTemp &ps)
std::string get_module_version()
static bool get_is_setup(Model *m, ParticleIndex pi)
The standard decorator for manipulating molecular structures.
Ints get_index(const ParticlesTemp &particles, const Subset &subset, const Subsets &excluded)
A decorator for a particle representing an atom.
void transform(XYZ a, const algebra::Transformation3D &tr)
Apply a transformation to the particle.
static Bonded setup_particle(Model *m, ParticleIndex pi)
void load_frame(RMF::FileConstHandle file, RMF::FrameID frame)
Load the given RMF frame into the state of the linked objects.
A decorator for a particle with x,y,z coordinates.
static Scale setup_particle(Model *m, ParticleIndex pi)
static Colored setup_particle(Model *m, ParticleIndex pi, Color color)
A decorator for a particle that is part of a rigid body but not rigid.
std::ostream & show(Hierarchy h, std::ostream &out=std::cout)
Print the hierarchy using a given decorator to display each node.
Find all nearby pairs by testing all pairs.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
def deprecated_object
Python decorator to mark a class as deprecated.
A decorator for a residue.
DensityGrid get_grid(IMP::em::DensityMap *in)
Return a dense grid containing the voxels of the passed density map.
General purpose algebraic and geometric methods that are expected to be used by a wide variety of IMP...
static bool get_is_setup(const IMP::ParticleAdaptor &p)
static bool get_is_setup(Model *m, ParticleIndex p)
Check if the particle has the needed attributes for a cast to succeed.
The general base class for IMP exceptions.
Rotation3D get_identity_rotation_3d()
Return a rotation that does not do anything.
void link_hierarchies(RMF::FileConstHandle fh, const atom::Hierarchies &hs)
Class to handle individual particles of a Model object.
Bond get_bond(Bonded a, Bonded b)
Get the bond between two particles.
std::string get_data_path(std::string file_name)
Return the full path to one of this module's data files.
int get_copy_index(Hierarchy h)
Walk up the hierarchy to find the current copy index.
double get_resolution(Model *m, ParticleIndex pi)
Estimate the resolution of the hierarchy as used by Representation.
A decorator for a rigid body.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Hierarchies get_leaves(const Selection &h)
A decorator for a molecule.
Select hierarchy particles identified by the biological name.
Support for the RMF file format for storing hierarchical molecular data and markup.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Transformation3D get_random_local_transformation(Vector3D origin, double max_translation=5., double max_angle_in_rad=0.26)
Get a local transformation.
Inferential scoring building on methods developed as part of the Inferential Structure Determination ...