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,-154182.69328281307,146.84014181906187
200,-150561.54412725402,202.21333959234298
300,-148167.39953218953,236.4575767860512
400,-146560.86349830756,258.436486685377
500,-145439.13934520682,272.4278095164674
600,-144398.90692638594,277.474536908782
700,-143922.56028575602,282.9090864511679
800,-143300.68101606035,284.0890902651323
900,-142839.77902464525,288.5596484705654
1000,-142368.53374238528,291.514018724144