Examples using CHARMM Files

There are popular ways to model biomolecular systems. Because solvation effects are often (always?) critically important to biological function, we need some way to model the solvent. The two popular approaches are to employ a continuum model with a dielectric constant equal to that of the bulk solvent or to model the solvent models directly in your system. These two approaches are termed _implicit_ and _explicit_ solvation, respectively.

The next 2 sections present examples using a Generalized Born implicit solvent model and explicit solvent based on the TIP3P water model. All of the files and examples here are included in the examples/charmm directory of the ParmEd release.

Generalized Born Implicit Solvent

For the purposes of this example, we are using an alanine pentapeptide. You can find the following files that you will need for this demonstration in the examples/charmm directory of the ParmEd distribution:

  • ala5_autopsf.psf

  • ala5_autopsf.pdb

The following sample script (simulate_charmm_gb.py in the ParmEd distribution) will set up and run the simulation using OpenMM:

#!/usr/bin/env python
from __future__ import division, print_function

import sys

# OpenMM Imports
import simtk.openmm as mm
import simtk.openmm.app as app

# ParmEd Imports
from parmed.charmm import CharmmPsfFile, CharmmCrdFile, CharmmParameterSet
from parmed.openmm import StateDataReporter
from parmed import unit as u

# Load the CHARMM files
print('Loading CHARMM files...')
params = CharmmParameterSet('toppar/par_all36_prot.prm')
ala5_gas = CharmmPsfFile('ala5_autopsf.psf')
ala5_crds = app.PDBFile('ala5_autopsf.pdb')

# NOTE NOTE
# The parameter set we used here is the CHARMM 36 force field, but this is
# strictly an example. It is important that you use the most accurate (typically
# most up-to-date) force fields for your own simulation. See the CHARMM
# parameter web page for updates:
# http://mackerell.umaryland.edu/CHARMM_ff_params.html
# END NOTE

# Create the OpenMM system
print('Creating OpenMM System')
system = ala5_gas.createSystem(params, nonbondedMethod=app.NoCutoff,
                               constraints=app.HBonds, implicitSolvent=app.HCT,
                               implicitSolventSaltConc=0.1*u.moles/u.liter,
)

# Create the integrator to do Langevin dynamics
integrator = mm.LangevinIntegrator(
                        300*u.kelvin,       # Temperature of heat bath
                        1.0/u.picoseconds,  # Friction coefficient
                        2.0*u.femtoseconds, # Time step
)

# Define the platform to use; CUDA, OpenCL, CPU, or Reference. Or do not specify
# the platform to use the default (fastest) platform
platform = mm.Platform.getPlatformByName('CUDA')
prop = dict(CudaPrecision='mixed') # Use mixed single/double precision

# Create the Simulation object
sim = app.Simulation(ala5_gas.topology, system, integrator, platform, prop)

# Set the particle positions
sim.context.setPositions(ala5_crds.positions)

# Minimize the energy
print('Minimizing energy')
sim.minimizeEnergy(maxIterations=500)

# Set up the reporters to report energies and coordinates every 100 steps
sim.reporters.append(
        StateDataReporter(sys.stdout, 100, step=True, potentialEnergy=True,
                          kineticEnergy=True, temperature=True,
                          volume=True, density=True)
)
sim.reporters.append(app.DCDReporter('ala5_gb.dcd', 100))

# Run dynamics
print('Running dynamics')
sim.step(10000)

Now we’ll dissect the script to help you understand what is happening at each step. We will divide the script into the sections following the print statements that announce when each stage begins.

Loading CHARMM files

In this stage, we instantiate the CharmmParameterSet with the parameters we wish to use. In this example we are using the CHARMM 36 force field, but be sure you check the CHARMM parameter website to make sure you use the latest and greatest force fields for your own work.

We also define the CharmmPsfFile object that stores all information about atom connectivity. The coordinates from PSF files generated by VMD are typically given in PDB format, so we use the OpenMM PDBFile object.

Create the OpenMM system

This command creates an OpenMM System object from the information stored in ala5_gas. It contains multiple Force instances for the bonds, angles, periodic torsions, improper torsions, Urey-Bradley forces, CMAP potentials, and nonbonded (electrostatic and van der Waals) interactions. It is in this function that we define the potential parameters we want to use. In this example, we have chosen the default values for each parameter except the ones specified. In particular:

  • params is the parameter set we defined: this is necessary to apply parameters to your system.

  • nonbondedMethod=app.NoCutoff indicates we do not want to use a cutoff for nonbonded interactions

  • constraints=app.HBonds indicates that we want to constrain all bonds in which at least one atom is hydrogen

  • implicitSolvent=app.HCT indicates that we want to use the Hawkins Cramer Truhlar GB model described in Hawkins et al., Chem. Phys. Lett., 1995, 51, p. 19824-19839.

  • implicitSolventSaltConc=0.1*u.liters/u.mole indicates that we want to model a ca. 0.1 molar solution of monovalent ions using a Debye screening model.

See the API documentation for the parmed package for a full listing of available options. If there are any other force objects you want to define, they can be added to the system after this step (like, for instance, positional restraints to a reference structure).

Create the integrator to do Langevin Dynamics

In this stage we specify an integrator. Common choices are LangevinIntegrator (as we’ve chosen here) to do simulations in the NVT ensemble and VerletIntegrator that allows us to do simulations either at constant energy or temperature if using the AndersenThermostat. In this example, we’ve chosen the Langevin integrator with a target temperature of 300 K, a friction coefficient of 1/ps and a time step of 2 fs.

Define the platform

In this stage, we define the platform we want to use. In this example we have chosen the CUDA platform, but this may not be available on every machine since it only runs on NVidia GPU hardware. Other choices are OpenCL (which will run on a variety of GPUs including those made by AMD/ATI and CPUs), CPU (which is an optimized version that runs natively on CPUs), and Reference (often quite slow).

The properties can be set for each platform. In this case, we specified that we wanted to use a mixed precision model (a good compromise between speed and precision).

Create the Simulation object

This step creates a Simulation object that will be used to run the actual simulations. If we wanted OpenMM to simply pick the fastest platform for us (rather than specify one directly), we could omit the platform and prop arguments.

Set the particle positions

This stage is very important. In this step, we set the particle positions stored in the ala5_crds object to our object. If you omit this step, you can get strange results or other errors like segmentation violations. These particle positions have been parsed from the PDB file, but you could use a CharmmCrdFile or CharmmRstFile (restart file) instead.

Minimize the energy

This stage performs a basic energy minimization to relax particle positions. This particular invocation will perform at most 500 iterations.

Set up the reporters

This stage defines reporters that will “report” on the status of the simulation periodically throughout the simulation. The first is an StateDataReporter which will print out a summary of energies and temperatures every 100 steps. Unlike the StateDataReporter in OpenMM, this reporter prints values in the AKMA unit system (Angstrom, Kilocalorie per mole, and atomic mass units). This reporter directs the printout to standard output (the screen), sys.stdout can be replaced with a different file-like object or a file name.

The second reporter is a DCD trajectory reporter, which is written in the standard DCD trajectory format.

Running dynamics

This is the stage that actually runs the MD. In this case, we are running 10,000 steps of MD. The wiki page with “Common recipes” has information regarding running a long simulation in chunks.

Explicit solvent

For the purposes of this example, we are using an alanine dipeptide. You can find the following files that you will need for this demonstration in the examples/charmm directory of the ParmEd distribution:

  • ala2_charmmgui.psf

  • ala2_charmmgui.pdb

The following sample script (simulate_charmm_gb.py in the ParmEd distribution) will set up and run the simulation using OpenMM:

#!/usr/bin/env python
from __future__ import division, print_function

import sys

# OpenMM Imports
import simtk.openmm as mm
import simtk.openmm.app as app

# ParmEd Imports
from parmed.charmm import CharmmPsfFile, CharmmCrdFile, CharmmParameterSet
from parmed.amber.openmmreporters import StateDataReporter
from parmed import unit as u

# Load the CHARMM files
print('Loading CHARMM files...')
params = CharmmParameterSet('toppar/par_all36_prot.prm',
                            'toppar/toppar_water_ions.str')
ala2_solv = CharmmPsfFile('ala2_charmmgui.psf')
ala2_crds = CharmmCrdFile('ala2_charmmgui.crd')

# NOTE NOTE
# The parameter set we used here is the CHARMM 36 force field, but this is
# strictly an example. It is important that you use the most accurate (typically
# most up-to-date) force fields for your own simulation. See the CHARMM
# parameter web page for updates:
# http://mackerell.umaryland.edu/CHARMM_ff_params.html
# END NOTE

# Compute the box dimensions from the coordinates and set the box lengths (only
# orthorhombic boxes are currently supported in OpenMM)
coords = ala2_crds.positions
min_crds = [coords[0][0], coords[0][1], coords[0][2]]
max_crds = [coords[0][0], coords[0][1], coords[0][2]]

for coord in coords:
    min_crds[0] = min(min_crds[0], coord[0])
    min_crds[1] = min(min_crds[1], coord[1])
    min_crds[2] = min(min_crds[2], coord[2])
    max_crds[0] = max(max_crds[0], coord[0])
    max_crds[1] = max(max_crds[1], coord[1])
    max_crds[2] = max(max_crds[2], coord[2])

ala2_solv.setBox(max_crds[0]-min_crds[0],
                 max_crds[1]-min_crds[1],
                 max_crds[2]-min_crds[2],
)

# Create the OpenMM system
print('Creating OpenMM System')
system = ala2_solv.createSystem(params, nonbondedMethod=app.PME,
                                nonbondedCutoff=12.0*u.angstroms,
                                constraints=app.HBonds,
                                switchDistance=10.0*u.angstroms,
)

# Create the integrator to do Langevin dynamics
integrator = mm.LangevinIntegrator(
                        300*u.kelvin,       # Temperature of heat bath
                        1.0/u.picoseconds,  # Friction coefficient
                        2.0*u.femtoseconds, # Time step
)

# Define the platform to use; CUDA, OpenCL, CPU, or Reference. Or do not specify
# the platform to use the default (fastest) platform
platform = mm.Platform.getPlatformByName('CUDA')
prop = dict(CudaPrecision='mixed') # Use mixed single/double precision

# Create the Simulation object
sim = app.Simulation(ala2_solv.topology, system, integrator, platform, prop)

# Set the particle positions
sim.context.setPositions(ala2_crds.positions)

# Minimize the energy
print('Minimizing energy')
sim.minimizeEnergy(maxIterations=500)

# Set up the reporters to report energies and coordinates every 100 steps
sim.reporters.append(
        StateDataReporter(sys.stdout, 100, step=True, potentialEnergy=True,
                          kineticEnergy=True, temperature=True,
                          volume=True, density=True)
)
sim.reporters.append(app.DCDReporter('ala2_solv.dcd', 100))

# Run dynamics
print('Running dynamics')
sim.step(10000)

Now we’ll dissect the script to help you understand what is happening at each step. We will divide the script into the sections following the print statements that announce when each stage begins.

Loading CHARMM files

In this stage, we instantiate the CharmmParameterSet with the parameters we wish to use. In this example we are using the CHARMM 36 force field, but be sure you check the CHARMM parameter website to make sure you use the latest and greatest force fields for your own work. Note that we also need to load the stream file with the water and ion parameters. The file type detection (RTF, parameter, or stream) is done by file-name extension (or the beginning of the file name if the file names end in .inp).

We also define the CharmmPsfFile object that stores all information about atom connectivity. The coordinates from PSF files generated by CHARMM GUI are given in the CHARMM coordinate format, so we use the CharmmCrdFile class to parse it.

Compute the box dimensions

In this step, we compute the box dimensions by calculating the extent of the coordinates in each of the X-, Y-, and Z-directions. Then the call to ala2_solv.setBox sets the periodic box vectors. You can optionally put the three angles between the lattice vectors (alpha, beta, and gamma) after the box lengths. The default angles are 90 degrees, since OpenMM only supports orthorhombic boxes.

Create the OpenMM system

This command creates an OpenMM System object from the information stored in ala2_solv. It contains multiple Force instances for the bonds, angles, periodic torsions, improper torsions, Urey-Bradley forces, CMAP potentials, and nonbonded (electrostatic and van der Waals) interactions. It is in this function that we define the potential parameters we want to use. In this example, we have chosen the default values for each parameter except the ones specified. In particular:

  • params is the parameter set we defined: this is necessary to apply parameters to your system.

  • nonbondedMethod=app.PME indicates we want to include the full electrostatic interactions using the Particle Mesh Ewald method

  • nonbondedCutoff=12.0*u.angstroms indicates we want to truncate van der Waals interactions and the direct-space sum in PME to 12 angstroms. (This is equivalent to the _ctofnb_ parameter in CHARMM)

  • constraints=app.HBonds indicates that we want to constrain all bonds in which at least one atom is hydrogen

  • switchDistance=10.0*u.angstroms indicates that we want to turn the switching function on at 10 Angstroms (This is equivalent to the ctonnb parameter in CHARMM)

See the API documentation for the parmed package for a full listing of available options (such as using a switching function for van der Waals interactions with the switchDistance keyword). If there are any other force objects you want to define, they can be added to the system after this step (like, for instance, positional restraints to a reference structure).

Create the integrator to do Langevin Dynamics

In this stage we specify an integrator. Common choices are LangevinIntegrator (as we’ve chosen here) to do simulations in the NVT ensemble and VerletIntegrator that allows us to do simulations either at constant energy or temperature if using the AndersenThermostat. In this example, we’ve chosen the Langevin integrator with a target temperature of 300 K, a friction coefficient of 1/ps and a time step of 2 fs.

Define the platform

In this stage, we define the platform we want to use. In this example we have chosen the CUDA platform, but this may not be available on every machine since it only runs on NVidia GPU hardware. Other choices are OpenCL (which will run on a variety of GPUs including those made by AMD/ATI and CPUs), CPU (which is an optimized version that runs natively on CPUs), and Reference (often quite slow).

The properties can be set for each platform. In this case, we specified that we wanted to use a mixed precision model (a good compromise between speed and precision).

Create the Simulation object

This step creates a Simulation object that will be used to run the actual simulations. If we wanted OpenMM to simply pick the fastest platform for us (rather than specify one directly), we could omit the platform and prop arguments.

Set the particle positions

This stage is very important. In this step, we set the particle positions stored in the ala2_crds object to our object. If you omit this step, you can get strange results or other errors like segmentation violations. These particle positions have been parsed from the CHARMM coordinate file, but you could use a PDBFile or CharmmRstFile (restart file) instead.

Minimize the energy

This stage performs a basic energy minimization to relax particle positions. This particular invocation will perform at most 500 iterations.

Set up the reporters

This stage defines reporters that will “report” on the status of the simulation periodically throughout the simulation. The first is an StateDataReporter which will print out a summary of energies and temperatures every 100 steps. Unlike the StateDataReporter in OpenMM, this reporter prints values in the AKMA unit system (Angstrom, Kilocalorie per mole, and atomic mass units). This reporter directs the printout to standard output (the screen), sys.stdout can be replaced with a different file-like object or a file name.

The second reporter is a DCD trajectory reporter, which is written in the standard DCD trajectory format.

Running dynamics

This is the stage that actually runs the MD. In this case, we are running 10,000 steps of MD. The wiki page with “Common recipes” has information regarding running a long simulation in chunks.