Skip to main content
The MD task runs atomistic molecular dynamics (MD) using ASE dynamics integrators. It supports:
  • Three statistical ensembles: NVE, NVT, NPT
  • Multiple thermostat/barostat algorithms
  • Time-dependent temperature and pressure schedules (annealing, ramps, etc.)
  • Trajectory checkpointing and restart
  • Dispersion corrections via TorchDFTD3Calculator
The implementation is adapted from the Atomate2 MLFF MD workflow.

Supported ensembles and dynamics

EnsembleDynamics stringASE class
nvevelocityverletase.md.verlet.VelocityVerlet
nvtlangevinase.md.langevin.Langevin
nvtandersenase.md.andersen.Andersen
nvtberendsenase.md.nvtberendsen.NVTBerendsen
nvtnose-hooverase.md.npt.NPT
nptnose-hooverase.md.npt.NPT
nptberendsenase.md.nptberendsen.NPTBerendsen
You can also pass any ASE MolecularDynamics class directly to dynamics.

Function signature

from mlip_arena.tasks import MD

result = MD(
    atoms=atoms,
    calculator=calc,
    ensemble="nvt",
    dynamics="langevin",
    time_step=None,
    total_time=1000,
    temperature=300.0,
    pressure=None,
    dynamics_kwargs=None,
    velocity_seed=None,
    zero_linear_momentum=True,
    zero_angular_momentum=True,
    traj_file=None,
    traj_interval=1,
    restart=True,
)

Parameters

body.atoms
ase.Atoms
required
Atomic structure to simulate. A copy is made internally.
body.calculator
ase.calculators.calculator.BaseCalculator
required
Calculator for energy and force evaluations at each MD step.
body.ensemble
string
default:"nvt"
Statistical ensemble. One of "nve", "nvt", or "npt".
body.dynamics
str | MolecularDynamics
default:"langevin"
Dynamics integrator. Accepts a string (see table above) or an ASE MolecularDynamics class. The string must be valid for the chosen ensemble.
body.time_step
number | None
default:"None"
Integration time step in femtoseconds. Defaults to 0.5 fs if hydrogen is present, otherwise 2.0 fs.
body.total_time
number
default:"1000"
Total simulation time in femtoseconds. The number of steps is int(total_time / time_step).
body.temperature
float | Sequence | ndarray | None
default:"300.0"
Target temperature in Kelvin. Can be:
  • A scalar for a constant temperature.
  • A sequence of values for a temperature schedule (linearly interpolated over the simulation).
  • Ignored for NVE.
body.pressure
float | Sequence | ndarray | None
default:"None"
Target pressure in eV/ų for NPT simulations. Can be a scalar or a sequence. Ignored for NVE and NVT.
body.dynamics_kwargs
dict | None
default:"None"
Extra keyword arguments forwarded to the ASE dynamics constructor. For Langevin, the default friction is 10.0 × 10⁻³ / fs (10 ps⁻¹) if not provided.
body.velocity_seed
integer | None
default:"None"
Random seed for Maxwell-Boltzmann velocity initialization. Set for reproducible simulations.
body.zero_linear_momentum
boolean
default:"true"
Remove the total linear momentum from the initial velocity distribution.
body.zero_angular_momentum
boolean
default:"true"
Remove the total angular momentum from the initial velocity distribution.
body.traj_file
str | Path | None
default:"None"
Path to an ASE .traj file for saving simulation frames. Parent directories are created automatically.
body.traj_interval
integer
default:"1"
Write a frame to traj_file every this many steps.
body.restart
boolean
default:"true"
If True and traj_file exists, resume the simulation from the last frame. Velocity and position information are restored from the trajectory.

Return value

KeyTypeDescription
atomsase.AtomsFinal structure after the simulation
runtimedatetime.timedeltaWall-clock time of the simulation
n_stepsintNumber of MD steps executed

Examples

from ase.build import bulk
from mlip_arena.tasks import MD
from mlip_arena.tasks.utils import get_calculator
from mlip_arena.models import MLIPEnum

atoms = bulk("Cu", "fcc", a=3.6) * (5, 5, 5)
calc = get_calculator(MLIPEnum.MACE_MP)

result = MD(
    atoms=atoms,
    calculator=calc,
    ensemble="nve",
    dynamics="velocityverlet",
    total_time=1000,   # 1 ps
    time_step=2.0,     # fs
    velocity_seed=42,
)

print(f"Steps: {result['n_steps']}")
print(f"Runtime: {result['runtime']}")

Temperature and pressure scheduling

Pass a list of values to temperature (or pressure for NPT) to define a piecewise schedule. The values are linearly interpolated over the total number of steps.
# Slow anneal from 300 K → 800 K → 300 K
temperature = [300, 800, 300]
The thermostat set-point is updated at every step via a callback, so the actual temperature tracks the schedule continuously.

Trajectory checkpointing

When traj_file is set and restart=True, the task detects an existing trajectory and continues from the last saved frame. This allows long simulations to be interrupted and resumed without losing progress.
For NPT simulations, ASE requires an upper-triangular cell. The task automatically applies a Schur decomposition to transform the cell before initializing the dynamics.