openmmml.models.deepmdpotential.DeepmdPotentialImpl#

class openmmml.models.deepmdpotential.DeepmdPotentialImpl(name, model: str = None, coordinatesCoefficient: float = 10.0, forceCoefficient: float = 964.8792534459, energyCoefficient: float = 96.48792534459)#

Implementation of the DeePMD (Deep Potential Molecular Dynamics) potential for OpenMM.

This class provides an interface to use DeePMD neural network potentials in OpenMM simulations through the OpenMMDeepmdPlugin. DeePMD potentials are trained using the DeepMD-kit framework and can provide highly accurate quantum mechanical-level forces and energies for molecular dynamics simulations.

The DeePMD potential supports various advanced features including: - Alchemical free energy calculations with lambda parameters - NNP/MM-style region selection - Custom unit conversions between different simulation packages

Prerequisites:
  • OpenMMDeepmdPlugin must be installed: conda install -c conda-forge ye-ding::openmm_deepmd_plugin

  • A trained DeePMD model file (.pb format from DeepMD-kit)

Usage:

The DeepmdPotentialImpl is typically used through the MLPotential interface:

from openmmml import MLPotential
import openmm.app as app

# Load your system
pdb = app.PDBFile('system.pdb')

# Create DeePMD potential
potential = MLPotential('deepmd', model='model.pb',
                        coordinatesCoefficient=10.0,
                        forceCoefficient=964.8792534459,
                        energyCoefficient=96.48792534459)

# Create system with DeePMD forces
system = potential.createSystem(
    topology=pdb.topology,
    atoms=None,  # Use all atoms, or specify subset
)

For alchemical simulations:

# Create system with lambda parameter for free energy calculations
potential = MLPotential('deepmd', model='model.pb')
system = potential.createSystem(
    topology=pdb.topology,
    lambdaName='lambda_alchemical',
    lambdaValue=0.5,  # Lambda value between 0 and 1
)

# During simulation, you can also update lambda value through the context:
context.setParameter('lambda_alchemical', new_lambda_value)

For NNP/MM-style simulations with atom subsets:

forcefield = ForceField('amber14-all.xml', 'amber14/tip3pfb.xml')
mm_system = forcefield.createSystem(topology)

# Define NNP region atoms (0-indexed)
nnp_atoms = [0, 1, 2, 3, 4]  # First 5 atoms

potential = MLPotential('deepmd', model='model.pb')
nnp_mm_system = potential.createMixedSystem(
    topology=pdb.topology,
    system=mm_system,
    atoms=nnp_atoms,  # Only these atoms use DeePMD
)

Note

The default unit conversion coefficients are set up for converting between DeePMD-kit’s native units (Angstrom, eV) and OpenMM’s units (nm, kJ/mol). Modify these coefficients if your model uses different units.

__init__(name, model: str = None, coordinatesCoefficient: float = 10.0, forceCoefficient: float = 964.8792534459, energyCoefficient: float = 96.48792534459) None#

Methods

__init__(name[, model, ...])

addForces(topology, system, atoms[, ...])

Add Force objects to a System to implement the potential function.

addForces(topology: Topology, system: System, atoms: Iterable[int] | None, forceGroup: int = 0, lambdaName: str | None = None, lambdaValue: float | None = 1.0, **args) None#

Add Force objects to a System to implement the potential function.

This is invoked by MLPotential.createSystem(). Subclasses must implement it to create the requested potential function.

Parameters:
  • topology (Topology) – the Topology from which the System is being created

  • system (System) – the System that is being created

  • atoms (Optional[Iterable[int]]) – the indices of atoms the potential should be applied to, or None if it should be applied to the entire System

  • forceGroup (int) – the force group that any newly added Forces should be in

  • args – any additional keyword arguments that were provided to createSystem() are passed to this method. This allows subclasses to customize their behavior based on extra arguments.