OpenMM
High-performance molecular dynamics simulation toolkit
Overview
OpenMM is a high-performance toolkit for molecular dynamics simulations, providing cross-platform support for CPU, GPU (CUDA, OpenCL), and multi-GPU acceleration. OpenMM is widely used in computational biology and chemistry for biomolecular simulations, drug discovery, and research applications.
Note
We provide OpenMM installation with optimised CPU platform support and comprehensive Python bindings. The installation is built with GCC 15.2.0 and includes full test suite validation.
Available versions
To view available openmm versions:
module avail openmm
Build recipes and configuration details are maintained in our GitLab repository:
Build optimizations
Our OpenMM installations are optimised for maximum performance on Discoverer’s hardware. We use GCC 15.2.0 to build OpenMM, as LLVM compilers have compatibility issues with OpenMM code.
Compiler optimizations:
RelWithDebInfo Build: Optimised code with debug symbols for performance and debugging support.
CPU-Specific Optimizations: Native architecture optimisations enabled, optimised for modern CPU architectures, Position Independent Code (PIC) for shared library support
Build system:
CMake: Modern build system configuration
Ninja: High-performance build tool for parallel compilation
SWIG: Python bindings generation
Parallel Compilation: Builds use parallel compilation jobs for efficient resource utilisation
Thread configuration:
Dynamic Thread Management: Thread count automatically adjusts based on available CPU cores
Conservative Limits: Thread limits prevent overparallelization during testing
CPU Affinity: CPU affinity binding enabled for optimal cache utilisation
OpenMP Integration: Multi-threaded CPU platform support
Testing:
Comprehensive Test Suite: All 133 tests pass successfully
Thread Configuration Validation: Thread limits verified during testing
Platform Functionality: CPU platform functionality verified
Warning
OpenMM code cannot handle Link Time Optimisation (LTO). LTO is disabled in our builds as it causes test failures. Additionally, LLVM compilers (versions 20 and 21) have bugs that cause illegal instruction exceptions, so we use GCC 15.2.0 instead.
Compiler support
Warning
OpenMM must be built with GCC 15.2.0 or later. LLVM compilers (versions 20 and 21) have bugs that cause illegal instruction exceptions in OpenMM code. LTO is not supported due to OpenMM code incompatibility.
Supported builds
Production builds:
module avail openmm
Available libraries
OpenMM provides shared libraries and Python bindings:
libOpenMM.so- OpenMM core libraryCore molecular dynamics simulation library providing force calculations, integrators, and platform abstraction.
Header files:
openmm/*.hLink flag:
-lOpenMMpkg-config:
OpenMM
libOpenMMCPU.so- CPU platform libraryCPU platform implementation with multi-threaded OpenMP support for molecular dynamics simulations.
Link flag:
-lOpenMMCPUPlatform: CPU (multi-threaded)
libOpenMMCuda.so- CUDA platform library (if CUDA available)CUDA platform implementation for GPU acceleration using NVIDIA GPUs.
Link flag:
-lOpenMMCudaPlatform: CUDA (GPU)
libOpenMMOpenCL.so- OpenCL platform library (if OpenCL available)OpenCL platform implementation for GPU acceleration using OpenCL-compatible devices.
Link flag:
-lOpenMMOpenCLPlatform: OpenCL (GPU)
Note
The libraries use optimised implementations and can be used in both C and C++ applications. Python bindings provide a user-friendly interface for most use cases.
Python integration
OpenMM provides comprehensive Python bindings that are automatically configured when you load the module. The Python bindings provide full access to all OpenMM functionality.
Loading the module:
# Load OpenMM module
module load openmm/8/<version>
# Verify Python bindings are available
python3 -c "import openmm; print(openmm.__version__)"
Basic usage example:
#!/usr/bin/env python3
import openmm
from openmm import app, unit
from openmm.app import PDBFile, ForceField, Simulation
from openmm import LangevinMiddleIntegrator, Platform
# Load PDB file
pdb = PDBFile('input.pdb')
# Create force field
forcefield = ForceField('amber14-all.xml', 'amber14/tip3pfb.xml')
# Create system
system = forcefield.createSystem(pdb.topology,
nonbondedMethod=app.PME,
nonbondedCutoff=1.0*unit.nanometer)
# Create integrator
integrator = LangevinMiddleIntegrator(300*unit.kelvin,
1.0/unit.picosecond,
0.002*unit.picosecond)
# Create simulation
simulation = Simulation(pdb.topology, system, integrator)
simulation.context.setPositions(pdb.positions)
# Minimize energy
simulation.minimizeEnergy()
# Run simulation
simulation.step(1000)
# Save final state
positions = simulation.context.getState(getPositions=True).getPositions()
app.PDBFile.writeFile(simulation.topology, positions, open('output.pdb', 'w'))
Platform selection:
import openmm
from openmm import Platform
# List available platforms
print("Available platforms:")
for i in range(Platform.getNumPlatforms()):
platform = Platform.getPlatform(i)
print(f" {i}: {platform.getName()}")
# Select CPU platform explicitly
platform = Platform.getPlatformByName('CPU')
# Create simulation with specific platform
simulation = Simulation(pdb.topology, system, integrator, platform)
Thread configuration:
import os
import openmm
# Set number of threads for CPU platform
os.environ['OPENMM_CPU_THREADS'] = '4'
os.environ['OMP_NUM_THREADS'] = '4'
os.environ['OMP_PROC_BIND'] = 'true'
os.environ['OMP_PLACES'] = 'cores'
# Create simulation (will use configured thread count)
simulation = Simulation(pdb.topology, system, integrator)
Example: Protein simulation with energy minimization:
#!/usr/bin/env python3
import openmm
from openmm import app, unit
from openmm.app import PDBFile, ForceField, Simulation
from openmm import LangevinMiddleIntegrator
# Load structure
pdb = PDBFile('protein.pdb')
# Create force field
forcefield = ForceField('amber14-all.xml', 'amber14/tip3pfb.xml')
# Create system with periodic boundary conditions
system = forcefield.createSystem(
pdb.topology,
nonbondedMethod=app.PME,
nonbondedCutoff=1.0*unit.nanometer,
constraints=app.HBonds
)
# Create integrator
integrator = LangevinMiddleIntegrator(
300*unit.kelvin, # Temperature
1.0/unit.picosecond, # Friction
0.002*unit.picosecond # Time step
)
# Create simulation
simulation = Simulation(pdb.topology, system, integrator)
simulation.context.setPositions(pdb.positions)
# Minimize energy
print("Minimizing energy...")
simulation.minimizeEnergy(maxIterations=1000)
# Run equilibration
print("Equilibrating...")
simulation.step(10000)
# Run production simulation
print("Running production simulation...")
simulation.step(100000)
# Save trajectory
positions = simulation.context.getState(getPositions=True).getPositions()
app.PDBFile.writeFile(simulation.topology, positions, open('output.pdb', 'w'))
# Get energy
state = simulation.context.getState(getEnergy=True)
print(f"Potential energy: {state.getPotentialEnergy()}")
Using OpenMM in SLURM jobs
OpenMM simulations are computationally intensive and should be run on compute nodes using SLURM batch jobs. The CPU platform uses OpenMP for multi-threading, so proper thread configuration and CPU binding are essential for optimal performance.
Basic SLURM job script:
#!/bin/bash
#SBATCH --partition=cn # On Discoverer CPU
#SBATCH --job-name=openmm_sim
#SBATCH --nodes=1
#SBATCH --ntasks-per-node=1
#SBATCH --cpus-per-task=24
#SBATCH --time=24:00:00
#SBATCH --mem=64G
#SBATCH --account=<your_project_slurm_account>
#SBATCH --qos=<your_project_slurm_qos>
#SBATCH --output=openmm_sim_%j.out
#SBATCH --error=openmm_sim_%j.err
module load openmm/8/<version>
module load python/<version>
# Run simulation
srun python3 simulation.py
CPU binding for optimal performance:
CPU binding is essential for OpenMM CPU platform simulations to prevent thread migration and improve cache locality:
#!/bin/bash
#SBATCH --partition=cn # On Discoverer CPU
#SBATCH --job-name=openmm_sim_cpu
#SBATCH --ntasks=1
#SBATCH --ntasks-per-node=1
#SBATCH --cpus-per-task=24
#SBATCH --time=24:00:00
#SBATCH --mem=64G
#SBATCH --account=<your_project_slurm_account>
#SBATCH --qos=<your_project_slurm_qos>
#SBATCH --output=openmm_cpu_%j.out
#SBATCH --error=openmm_cpu_%j.err
module load openmm/8/<version>
module load python/<version>
# Configure threads (use all allocated CPUs)
export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK
export OPENMM_CPU_THREADS=$SLURM_CPUS_PER_TASK
export OMP_PROC_BIND=true
export OMP_PLACES=cores
# Run with CPU binding
srun --cpu-bind=cores python3 run_simulation.py
Example: Production MD simulation:
#!/bin/bash
#SBATCH --partition=cn # On Discoverer CPU
#SBATCH --job-name=md_production
#SBATCH --ntasks=1
#SBATCH --cpus-per-task=24
#SBATCH --time=48:00:00
#SBATCH --mem=64G
#SBATCH --account=<your_project_slurm_account>
#SBATCH --qos=<your_project_slurm_qos>
#SBATCH --output=md_%j.out
#SBATCH --error=md_%j.err
module load openmm/8/<version>
module load python/<version>
# Thread configuration
export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK
export OPENMM_CPU_THREADS=$SLURM_CPUS_PER_TASK
export OMP_PROC_BIND=true
export OMP_PLACES=cores
# Run production simulation
srun --cpu-bind=cores python3 production_md.py
Example: Multiple independent simulations (array job):
#!/bin/bash
#SBATCH --partition=cn # On Discoverer CPU
#SBATCH --job-name=openmm_array
#SBATCH --array=1-10
#SBATCH --ntasks=1
#SBATCH --cpus-per-task=8
#SBATCH --time=12:00:00
#SBATCH --mem=16G
#SBATCH --account=<your_project_slurm_account>
#SBATCH --qos=<your_project_slurm_qos>
#SBATCH --output=sim_%A_%a.out
#SBATCH --error=sim_%A_%a.err
module load openmm/8/<version>
module load python/<version>
# Thread configuration
export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK
export OPENMM_CPU_THREADS=$SLURM_CPUS_PER_TASK
export OMP_PROC_BIND=true
export OMP_PLACES=cores
# Run simulation for array task
srun --cpu-bind=cores python3 run_simulation.py --input input_${SLURM_ARRAY_TASK_ID}.pdb
Example: GPU simulation (CUDA platform):
#!/bin/bash
#SBATCH --partition=cn # On Discoverer GPU
#SBATCH --job-name=openmm_gpu
#SBATCH --ntasks=1
#SBATCH --cpus-per-task=4
#SBATCH --gres=gpu:1
#SBATCH --time=24:00:00
#SBATCH --mem=32G
#SBATCH --account=<your_project_slurm_account>
#SBATCH --qos=<your_project_slurm_qos>
#SBATCH --output=openmm_gpu_%j.out
#SBATCH --error=openmm_gpu_%j.err
module load openmm/8/<version>
module load python/<version>
module load cuda/<version>
# Thread configuration (for CPU fallback if needed)
export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK
export OPENMM_CPU_THREADS=$SLURM_CPUS_PER_TASK
# Run GPU simulation
srun python3 gpu_simulation.py
Complete example: MD simulation script:
#!/usr/bin/env python3
"""
Example OpenMM molecular dynamics simulation script.
Usage: python3 md_simulation.py --input input.pdb --output output.pdb --steps 100000
"""
import argparse
import openmm
from openmm import app, unit
from openmm.app import PDBFile, ForceField, Simulation
from openmm import LangevinMiddleIntegrator, Platform
def main():
parser = argparse.ArgumentParser(description='Run MD simulation with OpenMM')
parser.add_argument('--input', required=True, help='Input PDB file')
parser.add_argument('--output', required=True, help='Output PDB file')
parser.add_argument('--steps', type=int, default=100000, help='Number of steps')
parser.add_argument('--temperature', type=float, default=300.0, help='Temperature (K)')
parser.add_argument('--platform', default='CPU', help='Platform (CPU, CUDA, OpenCL)')
args = parser.parse_args()
# Load structure
print(f"Loading structure from {args.input}")
pdb = PDBFile(args.input)
# Create force field
print("Creating force field...")
forcefield = ForceField('amber14-all.xml', 'amber14/tip3pfb.xml')
# Create system
print("Creating system...")
system = forcefield.createSystem(
pdb.topology,
nonbondedMethod=app.PME,
nonbondedCutoff=1.0*unit.nanometer,
constraints=app.HBonds
)
# Create integrator
integrator = LangevinMiddleIntegrator(
args.temperature*unit.kelvin,
1.0/unit.picosecond,
0.002*unit.picosecond
)
# Select platform
try:
platform = Platform.getPlatformByName(args.platform)
print(f"Using platform: {platform.getName()}")
except:
print(f"Platform {args.platform} not available, using CPU")
platform = Platform.getPlatformByName('CPU')
# Create simulation
simulation = Simulation(pdb.topology, system, integrator, platform)
simulation.context.setPositions(pdb.positions)
# Minimize energy
print("Minimizing energy...")
simulation.minimizeEnergy(maxIterations=1000)
# Run simulation
print(f"Running {args.steps} steps...")
simulation.step(args.steps)
# Save output
print(f"Saving output to {args.output}")
positions = simulation.context.getState(getPositions=True).getPositions()
app.PDBFile.writeFile(simulation.topology, positions, open(args.output, 'w'))
# Report energy
state = simulation.context.getState(getEnergy=True)
print(f"Final potential energy: {state.getPotentialEnergy()}")
print("Simulation complete!")
if __name__ == '__main__':
main()
SLURM job script for the example:
#!/bin/bash
#SBATCH --partition=cn # On Discoverer CPU
#SBATCH --job-name=md_example
#SBATCH --ntasks-per-node=1
#SBATCH --cpus-per-task=24
#SBATCH --time=24:00:00
#SBATCH --mem=64G
#SBATCH --account=<your_project_slurm_account>
#SBATCH --qos=<your_project_slurm_qos>
#SBATCH --output=md_example_%j.out
#SBATCH --error=md_example_%j.err
module load openmm/8/<version>
module load python/<version>
# Thread configuration
export OMP_NUM_THREADS=$SLURM_CPUS_PER_TASK
export OPENMM_CPU_THREADS=$SLURM_CPUS_PER_TASK
export OMP_PROC_BIND=true
export OMP_PLACES=cores
# Run simulation
srun --cpu-bind=cores python3 md_simulation.py \
--input input.pdb \
--output output.pdb \
--steps 1000000 \
--temperature 300.0 \
--platform CPU
Platform support
OpenMM supports multiple computational platforms for molecular dynamics simulations:
CPU Platform:
Multi-threaded execution with OpenMP
Optimised for modern CPU architectures
Suitable for small to medium-sized systems
Cross-platform compatibility
Thread count controlled by
OMP_NUM_THREADSandOPENMM_CPU_THREADS
CUDA Platform:
GPU acceleration using NVIDIA GPUs
Significant speedup for large systems
Multi-GPU support for very large simulations
Requires CUDA toolkit and compatible GPU
Automatically selected if available
OpenCL Platform:
GPU acceleration using OpenCL-compatible devices
Cross-vendor GPU support
Suitable for systems without CUDA support
Requires OpenCL runtime
Automatically selected if available
Reference Platform:
High-precision reference implementation
Useful for validation and debugging
Slower but more accurate
Always available
Platform selection:
import openmm
from openmm import Platform
# List available platforms
for i in range(Platform.getNumPlatforms()):
platform = Platform.getPlatformByName(Platform.getPlatform(i).getName())
print(f"{i}: {platform.getName()}")
# Select platform explicitly
platform = Platform.getPlatformByName('CPU') # or 'CUDA', 'OpenCL', 'Reference'
Performance considerations
CPU Platform:
Use CPU binding (
--cpu-bind=cores) to prevent thread migrationSet
OMP_NUM_THREADSandOPENMM_CPU_THREADSto match allocated CPUsEnable CPU affinity (
OMP_PROC_BIND=true,OMP_PLACES=cores)Allocate sufficient memory (typically 2-4 GB per 1000 atoms)
GPU Platform:
Use GPU platform for systems with >10,000 atoms
Allocate GPU resources (
--gres=gpu:1)CPU threads still needed for data transfer and coordination
Monitor GPU memory usage
Best practices:
Start with energy minimization before production runs
Use appropriate time steps (typically 0.002 ps for all-atom simulations)
Monitor energy conservation and system stability
Save checkpoints regularly for long simulations
Use appropriate force fields for your system type
Linking your application
If you need to link against OpenMM libraries in C or C++ code:
Method 1: Using environment variables (recommended)
# Load the module first
module load openmm/8/<version>
# Link against OpenMM - C++ code
g++ -o myapp myapp.cpp $CXXFLAGS $LDFLAGS -lOpenMM -lOpenMMCPU
Method 2: Using pkg-config
# Load the module first
module load openmm/8/<version>
# Link against OpenMM - C++ code
g++ -o myapp myapp.cpp $(pkg-config --cflags --libs OpenMM)
Method 3: Manual linking
# Load the module first
module load openmm/8/<version>
# Link against OpenMM - C++ code
g++ -o myapp myapp.cpp -I$OPENMM_ROOT/include -L$OPENMM_ROOT/lib -lOpenMM -lOpenMMCPU
Note
The Environment Modules automatically set CXXFLAGS and LDFLAGS when you load the module. Using these variables is the recommended approach as they remain correct even if the module path changes.
Getting help
For additional assistance:
See the Getting help documentation
OpenMM Documentation: https://docs.openmm.org
OpenMM Examples: https://github.com/openmm/openmm/tree/main/examples