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,-156790.456252678,118.08591831659236
200,-152916.5038932572,181.88421684850954
300,-149668.14625479927,219.46006314441277
400,-147357.36879939443,244.68090714252062
500,-146055.43195699117,260.04171601150136
600,-144690.5685493849,267.93570950286613
700,-144072.41516183136,281.97649060761586
800,-143119.16832046764,285.2726531984502
900,-142992.69708173777,290.75805435832194
1000,-142171.99426800164,294.0652531756243