1 from __future__
import print_function
17 from .
import atomicDominoUtilities
22 def __init__(self, parameterFileName):
23 self.createdMdOptimizer = 0
24 self.fixedNonFlexibleAtoms = 0
25 self.readParameters(parameterFileName)
26 self.wroteNativeProtein = 0
28 def loadHelpers(self):
29 self.loadDockingHelpers()
32 def readParameters(self, parameterFileName):
33 self.parameters = atomicDominoUtilities.readParameters(
37 def loadDockingHelpers(self):
38 outputDir = self.getParam(
"output_directory")
39 mkdirProcess = subprocess.Popen(
42 self.getParam(
'output_directory')],
44 stderr=subprocess.PIPE)
45 output = mkdirProcess.communicate()
47 self.namesToParticles = atomicDominoUtilities.makeNamesToParticles(
55 def initializeParticleNameDict(self, leaves):
59 name = self.quickParticleName(leaf)
60 particleDict[name] = {}
69 def createModel(self):
71 pdbName = self.getParam(
"native_pdb_input_file")
85 topology = ff.create_topology(self.protein)
86 topology.apply_default_patches()
88 topology.add_atom_types(self.protein)
90 ff.add_radii(self.protein)
91 ff.add_well_depths(self.protein)
93 self.topology = topology
96 def getParam(self, paramName):
97 paramValue = self.parameters[paramName]
101 def quickRestraintName(self, r):
103 ps = r.get_input_particles()
106 if (self.isAtomParticle(p) == 0):
109 if ((self.quickParticleName(p)
in self.flexibleAtoms) == 1):
111 nameList.append(
"%s_%s" % (self.quickParticleName(p), isFlexible))
113 return "%s acting on %s" % (r.get_name(),
" , ".join(nameList))
115 def quickParticleName(self, particle):
116 return atomicDominoUtilities.quickParticleName(particle)
119 def logTime(self, label):
120 nextTime = time.asctime(time.localtime(time.time()))
121 print(
"%s:\t %s" % (label, nextTime))
122 nextTimeEntry = [label, nextTime]
123 self.times.append(nextTimeEntry)
124 self.rawTimes.append(time.time())
127 def outputTimes(self):
128 print(
"Time for each step:")
129 for timeEntry
in self.times:
132 print(
"%s:\t%s" % (label, value))
134 def getFlexibleAtoms(self):
135 return self.flexibleAtoms
138 def getChain(self, particle):
142 chainId = chain.get_id()
146 def getResidue(self, particle):
149 residueNumber = residue.get_index()
152 def isAtomParticle(self, p):
153 return atomicDominoUtilities.isAtomParticle(p)
157 def allFlexibleAtoms(self, r):
158 ps = r.get_input_particles()
160 if (self.isAtomParticle(p) == 0):
162 if ((self.quickParticleName(p)
in self.flexibleAtoms) == 0):
169 def skipNonBondedPair(self, pair):
172 if (self.getResidue(pair[0]) == self.getResidue(pair[1])):
176 if ((self.quickParticleName(pair[0])
in self.flexibleAtoms) == 0
and (self.quickParticleName(pair[1])
in self.flexibleAtoms) == 0):
184 peptideChain = self.getParam(
"peptide_chain")
185 self.flexibleAtoms = {}
187 if (self.getParam(
"flexible_protein") ==
"yes"):
190 particlePairs = self.getClosePairs()
192 for pair
in particlePairs:
195 p0Chain = self.getChain(p0)
196 p1Chain = self.getChain(p1)
197 if (p0Chain == peptideChain
or p1Chain == peptideChain):
198 self.flexibleAtoms[self.quickParticleName(p0)] = 1
199 self.flexibleAtoms[self.quickParticleName(p1)] = 1
205 if (self.getChain(leaf) == peptideChain):
206 self.flexibleAtoms[self.quickParticleName(leaf)] = 1
207 self.fixNonFlexibleAtoms()
211 def getFixedAtoms(self):
215 leafName = self.quickParticleName(leaf)
216 if (leafName
in self.flexibleAtoms) == 0:
217 fixedAtoms[leafName] = 1
223 def initializeRandomExisting(self):
228 if ((self.quickParticleName(leaf)
in self.flexibleAtoms) == 1):
229 particlesToVary.append(self.quickParticleName(leaf))
230 atomToCoordinates = {}
231 for pName
in particlesToVary:
232 randomPname = random.choice(particlesToVary)
233 randomP = self.namesToParticles[randomPname]
234 print(
"setting coordinates for particle %s to the ones currently set for particle %s" % (pName, randomPname))
235 atomToCoordinates[pName] =
IMP.core.XYZ(randomP).get_coordinates()
237 for pName
in particlesToVary:
238 particle = self.namesToParticles[pName]
240 xyz.set_coordinates(atomToCoordinates[pName])
246 def setInitialPositions(self):
247 initialPositionMethod = self.getParam(
"initial_position_method")
252 if (initialPositionMethod ==
"random"):
254 myRandom = random.random()
255 initialOffset = float(
256 self.getParam(
"initial_atom_offset")) + myRandom
257 for particleName
in self.flexibleAtoms.keys():
258 p = self.namesToParticles[particleName]
260 coordinates = xyzDecorator.get_coordinates()
266 xyzDecorator.set_coordinates(randomXyz)
270 elif (initialPositionMethod ==
"random_existing"):
271 print(
"start random existing")
272 xyzs = self.initializeRandomExisting()
276 elif (initialPositionMethod ==
"random_existing_box"):
277 xyzs = self.initializeRandomExisting()
278 initialOffset = float(self.getParam(
"initial_atom_offset"))
279 for xyzDecorator
in xyzs:
280 coordinates = xyzDecorator.get_coordinates()
285 xyzDecorator.set_coordinates(randomXyz)
287 elif (initialPositionMethod ==
"random_box"):
295 for firstPname
in self.nonBondedPairs.keys():
296 secondPnames = self.nonBondedPairs[firstPname]
297 for secondPname
in secondPnames.keys():
298 if ((secondPname
in self.flexibleAtoms) == 0):
299 fixedAtom = self.namesToParticles[firstPname]
302 print(
"added next close atom %s to bounding box particles" % firstPname)
305 peptideParticles = []
308 if (self.getChain(leaf) == self.getParam(
"peptide_chain")):
313 fixedAtoms = self.getFixedAtoms()
314 for flexPname
in self.flexibleAtoms.keys():
316 flexParticle = self.namesToParticles[flexPname]
318 print(
"processing position for pname %s" % flexPname)
319 while(goodPosition == 0):
321 flexXyzDecorator.set_coordinates(
323 for fixedPname
in fixedAtoms.keys():
324 fixedParticle = self.namesToParticles[fixedPname]
329 print(
"finding new position for %s" % flexPname)
337 elif (initialPositionMethod ==
"random_full"):
339 for particleName
in self.flexibleAtoms.keys():
340 print(
"creating random position for particle %s" % particleName)
341 p = self.namesToParticles[particleName]
344 myRandom = random.random()
346 xyzDecorator.set_coordinates(randomXyz)
350 elif (initialPositionMethod ==
"file"):
352 initialPositionFile = self.getParam(
"saved_initial_atom_positions")
353 print(
"reading initial positions from %s" % initialPositionFile)
360 for initialLeaf
in initialLeaves:
361 name = self.quickParticleName(initialLeaf)
362 if ((name
in self.namesToParticles) == 1):
363 existingLeaf = self.namesToParticles[name]
366 existingLeafXyz.set_coordinates(
367 initialLeafXyz.get_coordinates())
369 print(
"Read in initial positions from file %s but this file contained particle %s which is not in current model" % (initialPositionFile, name))
371 print(
"Please specify valid initial position method (random or file)\n")
375 outputDir = self.getParam(
"output_directory")
376 initialAtomPositionFile = self.getParam(
377 "initial_atom_positions_output_file")
378 fullOutputFile = os.path.join(outputDir, initialAtomPositionFile)
379 print(
"writing output to file %s" % fullOutputFile)
383 initialScore = self.model.evaluate(
False)
384 print(
"initial score for model is %s" % initialScore)
387 def getNonBondedPairs(self):
388 nonBondedDefinition = self.getParam(
"non_bonded_definition")
390 if (nonBondedDefinition ==
"close_pairs"):
391 particlePairs = self.getClosePairs()
392 elif (nonBondedDefinition ==
"random"):
393 particlePairs = self.getRandomParticlePairs()
394 elif (nonBondedDefinition ==
"manual"):
395 particlePairs = self.readManualRestraints(nativeProtein, 6)
397 print(
"Please specify valid non bonded pair definition")
402 def getXyzDistance(self, pair):
410 def getClosePairs(self):
415 closePairDistance = float(self.getParam(
"close_pair_distance"))
418 useHardCutoff = self.getParam(
"use_hard_distance_cutoff")
422 self.model.evaluate(
False)
424 particlePairs = nbl.get_particle_pairs()
426 for pair
in particlePairs:
427 distance = self.getXyzDistance(pair)
428 q0 = self.quickParticleName(pair[0])
429 q1 = self.quickParticleName(pair[1])
431 if (useHardCutoff ==
"0" or (useHardCutoff ==
"1" and distance < closePairDistance)):
432 finalPairs.append([pair[0], pair[1]])
437 def addForceFieldRestraints(self):
439 useForcefieldRestraints = self.getParam(
"use_forcefield_restraints")
440 if (useForcefieldRestraints ==
"1"):
448 bonds = self.topology.add_bonds(self.protein)
449 angles = ff.create_angles(bonds)
450 dihedrals = ff.create_dihedrals(bonds)
451 impropers = self.topology.add_impropers(self.protein)
455 self.ffParticleGroups = self.initializeParticleNameDict(
458 self.createForceFieldRestraintsForModel(
466 self.createForceFieldRestraintsForModel(
474 self.createForceFieldRestraintsForModel(
477 "Dihedral Restraints",
479 self.createForceFieldRestraintsForModel(
482 "Improper Restraints",
489 print (
"Skipping forcefield restraints")
493 def createForceFieldRestraintsForModel(
500 maxForceFieldScore = float(self.getParam(
"max_forcefield_score"))
502 cont.add_particles(ffParticles)
504 singletonScore, cont, restraintName)
508 decomposedRestraintTemp = listBondRestraint.create_decomposition()
509 decomposedRestraints = IMP.RestraintSet.get_from(
510 decomposedRestraintTemp)
511 rs_restraints = decomposedRestraints.get_restraints()
514 del listBondRestraint
516 count = decomposedRestraints.get_number_of_restraints()
517 for i
in range(count):
519 r = decomposedRestraints.get_restraint(i)
521 if (self.allFlexibleAtoms(r) == 1):
523 self.model.add_restraint(r)
524 r.set_maximum_score(maxForceFieldScore)
525 self.addParticlesToGroups(r)
533 def addParticlesToGroups(self, r):
535 ps = r.get_input_particles()
539 if (self.isAtomParticle(p) == 0):
541 inputParticles.append(p.get_name())
543 for firstName
in inputParticles:
544 for secondName
in inputParticles:
545 if (firstName == secondName):
547 self.ffParticleGroups[firstName][secondName] = 1
552 def addExcludedVolumeRestraints(self):
554 restrainedParticles = atomicDominoUtilities.getRestrainedParticles(
555 self.protein, self.model, self.namesToParticles)
559 self.model.add_restraint(evr)
564 def addFlexFlexLjRestraints(self, scalingFactor):
565 if (self.initializedFlexFlexLjRestraints == 0):
566 print(
"adding initial lennard jones restraints between pairs of flexible atoms")
567 self.initializedFlexFlexLjRestraints = 1
573 for i
in range(len(leaves)):
575 for j
in range(i + 1, len(leaves)):
577 if ((self.quickParticleName(iLeaf)
in self.flexibleAtoms) == 1
and (self.quickParticleName(jLeaf)
in self.flexibleAtoms) == 1):
579 ps, [iLeaf, jLeaf],
"LennardJones %s_%s" %
581 self.model.add_restraint(ljpr)
583 self.protein, IMP.atom.CHAIN_TYPE)
584 peptideChainHierarchy =
None
585 print(
"scaling atomic radii by scaling factor %s" % scalingFactor)
587 self.ff.add_radii(self.protein)
588 for pName
in self.flexibleAtoms.keys():
589 particle = self.namesToParticles[pName]
591 radius = xyzr.get_radius()
592 radius *= scalingFactor
593 xyzr.set_radius(radius)
597 def addFixedFlexibleLjRestraints(self):
598 print(
"adding initial fixed flexible lj restraints")
599 fixedAtoms = self.getFixedAtoms()
603 for flexPname
in self.flexibleAtoms.keys():
604 flexParticle = self.namesToParticles[flexPname]
605 for fixedPname
in fixedAtoms.keys():
606 fixedParticle = self.namesToParticles[fixedPname]
608 ps, [fixedParticle, flexParticle],
"LennardJones %s_%s" %
609 (flexPname, fixedPname))
610 self.model.add_restraint(ljpr)
613 def createDopeRestraint(self, pair):
614 q0 = self.quickParticleName(pair[0])
615 q1 = self.quickParticleName(pair[1])
617 prName =
"NonBonded %s_%s" % (q0, q1)
619 self.dopePairScore, pair, prName)
622 self.model.add_restraint(pairRestraint)
623 maxNonBondedScore = float(self.getParam(
"max_non_bonded_score"))
624 pairRestraint.set_maximum_score(maxNonBondedScore)
629 def addClosePairNonBondedRestraints(self):
630 particlePairs = self.getNonBondedPairs()
632 nonBondedRestraintType = self.getParam(
"non_bonded_restraint_type")
633 dopeCutoff = float(self.getParam(
"dope_cutoff"))
638 if (nonBondedRestraintType ==
"dope"):
641 self.dopePairScore = ps
643 self.nonBondedPairs = self.initializeParticleNameDict(
645 for pair
in particlePairs:
647 if (self.skipNonBondedPair(pair) == 1):
650 self.createDopeRestraint(pair)
653 self.nonBondedPairs[self.quickParticleName(pair[0])][
654 self.quickParticleName(pair[0])] = 1
655 self.nonBondedPairs[self.quickParticleName(pair[1])][
656 self.quickParticleName(pair[1])] = 1
658 print(
"Added %s non-bonded restraints" % nonBondedCount)
659 self.model.evaluate(
False)
661 def addCompleteNonBondedRestraints(self):
665 for i
in range(len(leaves)):
667 for j
in range(i + 1, len(leaves)):
670 if ((self.quickParticleName(jLeaf)
in self.nonBondedPairs[self.quickParticleName(iLeaf)]) == 1):
672 if (self.skipNonBondedPair([iLeaf, jLeaf]) == 1):
675 self.createDopeRestraint([iLeaf, jLeaf])
680 print(
"added %s restraints for %s particles" % (counter, len(leaves)))
681 self.model.evaluate(
False)
684 def writeAllRestraints(self):
685 print(
"write all restraints:")
687 print(self.quickRestraintName(r))
690 def initializeMdSchedule(self):
691 schedule = self.getParam(
"md_schedule")
692 stages = schedule.split()
695 print(
"unpacking stage %s" % stage)
696 [steps, temperature, excludedVolume] = stage.split(
'_')
697 totalSteps += int(steps)
698 self.totalMdSteps = totalSteps
699 self.mdStages = stages
704 def createMdRootHandle(self):
707 outputDir = self.getParam(
"output_directory")
709 if (self.getParam(
"read_assignments_from_file") ==
"no"):
711 mdFileName = self.getParam(
"md_trajectory_output_file")
712 fileName = os.path.join(outputDir, mdFileName)
713 print(
"creating rmf file with filename %s" % fileName)
714 rh = RMF.create_rmf_file(fileName)
715 my_kc = rh.add_category(
"my data")
721 mdFileName = self.getParam(
"saved_md_filename")
723 rh = RMF.open_rmf_file(mdFileName)
724 IMP.rmf.set_hierarchies(rh, [self.protein])
728 def evaluateModel(self):
729 return self.model.evaluate(
False)
731 def loadMdHelpers(self):
733 self.createMdRootHandle()
735 self.createMdOptimizer()
737 self.addMdOptimizerStates()
739 self.initializeMdSchedule()
741 self.createdMdOptimizer = 1
744 def createMdOptimizer(self):
745 print(
"creating md object")
750 md.set_velocity_cap(100.0)
751 self.mdOptimizer = md
753 def addMdOptimizerStates(self):
759 vsIndex = self.mdOptimizer.add_optimizer_state(vsos)
762 hdos = IMP.rmf.SaveHierarchyConfigurationOptimizerState(
763 [self.protein], self.rootHandle)
764 self.mdOptimizer.add_optimizer_state(hdos)
769 def readTrajectories(self):
771 cgFileName = self.getParam(
"cg_output_file")
772 bestCgScore = 10000000
774 bestCgRmsd = 10000000
777 outputDir = self.getParam(
"output_directory")
778 trajectoryFile = self.getParam(
"md_trajectory_output_file")
779 fullFile = os.path.join(outputDir, trajectoryFile)
780 rh = RMF.open_rmf_file(fullFile)
781 IMP.rmf.set_hierarchies(rh, [self.protein])
782 framesToRead = atomicDominoUtilities.getMdIntervalFrames(
783 rh, int(self.getParam(
"cg_interval")), self.protein)
785 if (len(framesToRead) > 0):
786 for cgNumber
in framesToRead:
788 outputDir = self.getParam(
"output_directory")
789 fullCgFileName = os.path.join(
791 (cgFileName, cgNumber))
792 rh = RMF.open_rmf_file(fullCgFileName)
793 IMP.rmf.set_hierarchies(rh, [self.protein])
795 frameCount = IMP.rmf.get_number_of_frames(rh, self.protein)
797 score = self.model.evaluate(
False)
798 rmsd = self.calculateNativeRmsd(self.flexibleAtoms)
799 print(
"cg number %s has score %s rmsd %s" % (cgNumber, score, rmsd))
800 if (score < bestCgScore):
802 bestCgScoreFile = fullCgFileName
803 if (rmsd < bestCgRmsd):
805 bestCgRmsdFile = fullCgFileName
808 self.singlePdbResults(
811 self.getParam(
"best_cg_score_output_file"))
812 self.singlePdbResults(
815 self.getParam(
"best_cg_rmsd_output_file"))
816 self.singlePdbResults(
818 (cgFileName, framesToRead[-1]), -1, self.getParam(
"final_cg_frame_output_file"))
819 finalCgRmsd = self.calculateNativeRmsd(self.flexibleAtoms)
820 self.bestCgScore = bestCgScore
821 self.bestCgRmsd = bestCgRmsd
822 self.finalCgRmsd = finalCgRmsd
824 def singlePdbResults(self, trajectoryFile, frame, outputPdbFile):
826 fullTrajectoryFile = os.path.join(
827 self.getParam(
"output_directory"),
829 fullOutputFile = os.path.join(
830 self.getParam(
"output_directory"),
832 rh = RMF.open_rmf_file(fullTrajectoryFile)
833 IMP.rmf.set_hierarchies(rh, [self.protein])
835 frame = IMP.rmf.get_number_of_frames(rh, self.protein) - 1
839 def calculateRmsd(self, otherProtein, flexibleAtoms):
840 otherNamesToParticles = atomicDominoUtilities.makeNamesToParticles(
844 for pName
in otherNamesToParticles.keys():
845 if ((pName
in flexibleAtoms) == 0):
847 otherParticle = otherNamesToParticles[pName]
848 modelParticle = self.namesToParticles[pName]
856 def calculateNativeRmsd(self, flexibleAtoms):
858 if (self.wroteNativeProtein == 0):
859 pdbName = self.getParam(
"native_pdb_input_file")
865 self.wroteNativeProtein = 1
867 return self.calculateRmsd(self.nativeProtein, flexibleAtoms)
870 def runMolecularDynamics(self):
871 self.initializedFlexFlexLjRestraints = 0
872 if (self.fixedNonFlexibleAtoms == 0):
873 print(
"error: before running md, must fix non flexible atoms")
876 for stage
in self.mdStages:
877 [steps, temperature, excludedVolume] = stage.split(
'_')
878 stepsSoFar += int(steps)
879 self.vsos.set_temperature(int(temperature))
880 if (excludedVolume ==
"1"):
881 self.addFixedFlexibleLjRestraints()
882 elif (excludedVolume !=
"0"):
884 if (self.getParam(
"excluded_volume_type") ==
"lj"):
885 self.addFlexFlexLjRestraints(float(excludedVolume))
886 elif (self.getParam(
"excluded_volume_type") ==
"ev"):
887 self.addExcludedVolumeRestraints(
890 "please set excluded_volume_type to either lj or ev"
892 print(
"running md at temperature %s for %s steps" % (temperature, steps))
893 self.mdOptimizer.optimize(int(steps))
902 print(
"model score after stage %s is %s" % (stage, self.model.evaluate(
False)))
903 print(
"done running md")
906 cgInterval = int(self.getParam(
"cg_interval"))
907 outputDir = self.getParam(
"output_directory")
908 trajectoryFile = self.getParam(
"md_trajectory_output_file")
909 fullFile = os.path.join(outputDir, trajectoryFile)
910 print(
"open rmf %s" % fullFile)
911 rh = RMF.open_rmf_file(fullFile)
912 IMP.rmf.set_hierarchies(rh, [self.protein])
913 framesToRead = atomicDominoUtilities.getMdIntervalFrames(
914 rh, cgInterval, self.protein)
916 for i
in framesToRead:
917 cgFileName =
"%s%s" % (self.getParam(
"cg_output_file"), i)
918 self.applyCg(i, cgFileName)
922 def applyCg(self, mdStepStart, cgFileName):
924 cgSteps = int(self.getParam(
"cg_steps"))
927 print(
"apply cg: loading md frame for cg step start %s" % mdStepStart)
931 outputDir = self.getParam(
"output_directory")
932 fileName = os.path.join(outputDir, cgFileName)
933 print(
"creating cg hdf5 file %s" % fileName)
934 rh = RMF.create_rmf_file(fileName)
935 my_kc = rh.add_category(
"my data")
940 firstScore = self.model.evaluate(
False)
942 print(
"creating optimizer state")
943 hdos = IMP.rmf.SaveHierarchyConfigurationOptimizerState(
945 hdos.set_skip_steps(0)
947 cg.add_optimizer_state(hdos)
951 secondScore = self.model.evaluate(
False)
952 print(
"cg score after md step %s before %s after %s" % (mdStepStart, firstScore, secondScore))
956 def fixNonFlexibleAtoms(self):
957 if (self.createdMdOptimizer == 0):
958 self.createMdOptimizer()
962 self.fixedNonFlexibleAtoms = 1
965 if ((self.quickParticleName(leaf)
in self.flexibleAtoms) == 0):
967 xyzDecorator.set_coordinates_are_optimized(0)
969 def writeOutput(self):
971 flexiblePeptideAtoms = 0
972 flexibleProteinAtoms = 0
973 fixedPeptideAtoms = 0
974 fixedProteinAtoms = 0
977 if ((self.quickParticleName(leaf)
in self.flexibleAtoms) == 0):
978 if (self.getChain(leaf) == self.getParam(
"peptide_chain")):
979 fixedPeptideAtoms += 1
981 fixedProteinAtoms += 1
983 if (self.getChain(leaf) == self.getParam(
"peptide_chain")):
984 flexiblePeptideAtoms += 1
986 flexibleProteinAtoms += 1
987 print(
"Flexible peptide atoms: %s" % flexiblePeptideAtoms)
988 print(
"Fixed peptide atoms: %s" % fixedPeptideAtoms)
989 print(
"Flexible protein atoms: %s" % flexibleProteinAtoms)
990 print(
"Fixed protein atoms: %s" % fixedProteinAtoms)
992 def writeOsOutput(self):
994 print(
"Best cg score: %s" % self.bestCgScore)
995 print(
"best cg rmsd: %s" % self.bestCgRmsd)
996 print(
"final cg rmsd: %s" % self.finalCgRmsd)
Applies a SingletonScore to each Singleton in a list.
Chain get_chain(Hierarchy h)
Score the angle based on a UnaryFunction,.
Simple conjugate gradients optimizer.
algebra::BoundingBoxD< 3 > get_bounding_box(const Hierarchy &h)
Get a bounding box for the Hierarchy.
void write_pdb(const Selection &mhd, TextOutput out, unsigned int model=1)
Hierarchy get_residue(Hierarchy mhd, unsigned int index)
Get the residue with the specified index.
Return all close unordered pairs of particles taken from the SingletonContainer.
double get_distance(XYZR a, XYZR b)
Compute the sphere distance between a and b.
void read_pdb(TextInput input, int model, Hierarchy h)
Vector3D get_random_vector_in(const Cylinder3D &c)
Generate a random vector in a cylinder with uniform density.
CHARMM force field parameters.
static Float get_k_from_standard_deviation(Float sd, Float t=297.15)
Return the k to use for a given Gaussian standard deviation.
Class for storing model, its restraints, constraints, and particles.
Maintains temperature during molecular dynamics by velocity scaling.
Select all non-alternative ATOM records.
double get_rmsd(const Selection &s0, const Selection &s1)
void add_hierarchy(RMF::FileHandle fh, atom::Hierarchy hs)
Lennard-Jones score between a pair of particles.
Simple molecular dynamics optimizer.
Store a list of ParticleIndexes.
A decorator for a particle representing an atom.
Hierarchies get_by_type(Hierarchy mhd, GetByType t)
Gather all the molecular particles of a certain level in the hierarchy.
algebra::BoundingBoxD< 3 > get_bounding_box(const XYZRs &ps)
Get the bounding box.
Score the improper dihedral based on a UnaryFunction,.
RestraintsTemp get_restraints(It b, It e)
void load_frame(RMF::FileConstHandle file, RMF::FrameID frame)
A decorator for a particle with x,y,z coordinates.
Score the bond based on a UnaryFunction,.
Basic functionality that is expected to be used by a wide variety of IMP users.
General purpose algebraic and geometric methods that are expected to be used by a wide variety of IMP...
Prevent a set of particles and rigid bodies from inter-penetrating.
void add_dope_score_data(atom::Hierarchy h)
double get_distance(const VectorD< D > &v1, const VectorD< D > &v2)
Compute the distance between two vectors.
Applies a PairScore to a Pair.
std::string get_data_path(std::string file_name)
Return the full path to one of this module's data files.
Smooth interaction scores by switching the derivatives (force switch).
Hierarchies get_leaves(const Selection &h)
Support for the RMF file format for storing hierarchical molecular data and markup.
Divide-and-conquer inferential optimization in discrete space.
Score the dihedral angle.
Harmonic function (symmetric about the mean)
A decorator for a particle with x,y,z coordinates and a radius.