First SimulationΒΆ
This is a good starting point for learning how to use OpenMM.
It loads a PDB file villin.pdb, which is the villin headpiece protein in a box of water. In then parameterizes it using the Amber14 forcefield and TIP3P-FB water model, energy minimizes it, and simulates it for 1000 steps with a langevin intergator.
[1]:
from openmm.app import *
from openmm import *
from openmm.unit import *
from sys import stdout
# Load in the PDB strucure
pdb = PDBFile('villin.pdb')
# Specifiy the forcefield
forcefield = ForceField('amber14-all.xml', 'amber14/tip3pfb.xml')
# Combine the molecular topology and the forcefield
system = forcefield.createSystem(pdb.topology, nonbondedMethod=PME,
nonbondedCutoff=1*nanometer, constraints=HBonds)
# Create the integrator to use for advacing the equations of motion.
# It specifies a Langevin integrator.
# The paramters set are temperature, friction coefficient, and timestep.
integrator = LangevinMiddleIntegrator(300*kelvin, 1/picosecond, 0.004*picoseconds)
# Combines the molecular topology, system, and integrator
# to begin a new simulation.
simulation = Simulation(pdb.topology, system, integrator)
simulation.context.setPositions(pdb.positions)
# Perform local energy minimization
print("Minimizing energy...")
simulation.minimizeEnergy(maxIterations=100)
# Write the trajectory to a file called "output.pdb"
simulation.reporters.append(PDBReporter('output.pdb', 1000))
# Report infomation to the screen as the simulation runs
simulation.reporters.append(StateDataReporter(stdout, 100, step=True,
potentialEnergy=True, temperature=True))
#Run the simulation for 1000 timsteps
print("Running simulation...")
simulation.step(1000)
Minimizing energy...
Running simulation...
#"Step","Potential Energy (kJ/mole)","Temperature (K)"
100,-154068.74212126824,152.4682082389034
200,-150566.10974992323,203.50078962899158
300,-147603.0766119699,232.42649677995846
400,-145960.55765141113,252.44603316954334
500,-144946.2349824087,267.4488476759998
600,-143932.46107307085,279.74335768872004
700,-143353.40083165228,284.42324331455455
800,-142617.5256855465,287.16801921337355
900,-142302.41056646284,287.73885604690383
1000,-141892.96290554863,294.4363570473461