Molecular Dynamics
Overview
This module provides functionalities for running molecular dynamics (MD) simulations using ASE (Atomic Simulation Environment). It includes different thermostats and integrators to handle NVT MD simulations.
Each iteration of ML/MD data will be stored in corresponding iter_0000xx directory.
Which will have the a corresponding log_file and traj file.
>>> cat Iter1_dpmd.log
Time[ps] Etot[eV] Epot[eV] Ekin[eV] T[K]
0.0000 -112.0807 -112.8950 0.8143 300.0
0.0700 -111.6322 -112.7149 1.0828 398.9
0.1400 -112.4215 -113.3518 0.9303 342.7
0.2100 -112.9996 -113.6775 0.6779 249.8
0.2800 -112.6910 -113.7220 1.0310 379.8
0.3500 -112.8007 -113.2903 0.4896 180.4
Features:
Nose-Hoover and Langevin thermostats for NVT simulations
ab-initio energy calculations
ab-initio and ML molecular dynamics execution
Usage Examples
The examples below correspond to the aimd_setup (or mlip_setup)
blocks in input.yaml. The Python API mirrors those settings directly.
NVT — Nose-Hoover thermostat
aimd_setup:
ensemble: "NVT"
temperature: 300.0
timestep_fs: 1.0
steps: 500
thermostat:
type: "Nose"
tdamp: 2.0
from ase.io import read
from sparc.src.ase_md import NoseNVT
atoms = read("POSCAR")
dyn = NoseNVT(atoms, temperature=300, tdamp=2.0, timestep=1.0)
dyn.run(500)
NVT — Langevin thermostat
aimd_setup:
ensemble: "NVT"
temperature: 300.0
timestep_fs: 1.0
steps: 500
thermostat:
type: "Langevin"
friction: 0.01
from sparc.src.ase_md import LangevinNVT
dyn = LangevinNVT(atoms, temperature=300, friction=0.01, timestep=1.0)
dyn.run(500)
NVT — Temperature ramp (Langevin)
Temperature ramping linearly from temperature to temp_end over the
run. Nose-Hoover resists rapid temperature changes; use Langevin for ramping.
aimd_setup:
ensemble: "NVT"
temperature: 300.0
temp_end: 1000.0
timestep_fs: 1.0
steps: 1000
thermostat:
type: "Langevin"
friction: 0.01
from sparc.src.ase_md import LangevinNVT
dyn = LangevinNVT(atoms, temperature=300, temp_end=1000,
friction=0.01, timestep=1.0)
dyn.run(1000)
NVE — microcanonical ensemble
aimd_setup:
ensemble: "NVE"
temperature: 300.0
timestep_fs: 1.0
steps: 500
from sparc.src.ase_md import ExecuteAbInitioDynamics
from ase.md.verlet import VelocityVerlet
import ase.units as units
dyn = VelocityVerlet(atoms, timestep=1.0 * units.fs)
ExecuteAbInitioDynamics(system=atoms, dyn=dyn, steps=500, pace=10,
log_filename="aimd.log", trajfile="aimd.traj",
dir_name="aimd_results", name="NVE_AIMD")
NPT — isothermal-isobaric ensemble
Pressure is given in bar and compressibility in 1/bar, following
ASE’s units.bar convention (pressure_au = pressure * units.bar,
compressibility_au = compressibility / units.bar).
aimd_setup:
ensemble: "NPT"
temperature: 300.0
timestep_fs: 1.0
steps: 500
tau_t: 100.0 # Berendsen temperature coupling time [fs]
tau_p: 1000.0 # Berendsen pressure coupling time [fs]
pressure: 1.01325 # Target pressure [bar] (1 atm = 1.01325 bar)
compressibility: 4.57e-5 # Isothermal compressibility [1/bar] (optional; default: Cu ~7.1e-7 bar⁻¹)
Note
For NPT, the thermostat.type key is ignored. NPTBerendsen
handles both temperature and pressure coupling via tau_t and tau_p.
from sparc.src.ase_md import NPT
dyn = NPT(atoms, timestep=1.0, temperature=300,
tau_t=100.0, pressure=1.01325, tau_p=1000.0,
compressibility=4.57e-5)
Parameter |
Unit |
Notes |
|---|---|---|
|
bar |
1 atm ≈ 1.01325 bar; 1 GPa = 10 000 bar |
|
1/bar |
Default: 7.1×10⁻⁷ bar⁻¹ (Cu). Water ≈ 4.57×10⁻⁵ bar⁻¹ |
|
fs |
Berendsen temperature coupling time constant |
|
fs |
Berendsen pressure coupling time constant |
Module Contents
ASE-based Molecular Dynamics module for SPARC package.
This module provides MD simulation capabilities using ASE’s MD integrators, supporting various ensembles (NVE, NVT, NPT) with different thermostats.
- sparc.src.ase_md.CalculateDFTEnergy(idx, header, system, log_filename, dir_name, trajfile)[source]
Calculate the DFT energy and forces for a candidate structure.
- Parameters:
idx (
int) – Candidate index (used in log output)header (
bool) – If True, write column header to log filesystem (
Atoms) – Candidate structure with DFT calculator attachedlog_filename (
str) – Filename for the per-candidate energy logdir_name (
str) – Directory where log and trajectory files are writtentrajfile (
str) – Filename for the ASE trajectory file
- sparc.src.ase_md.ExecuteAbInitioDynamics(system, dyn, steps, pace, log_filename, trajfile, dir_name, name, temp_start=None, temp_end=None)[source]
Run an ab initio MD simulation with DFT calculator.
- Parameters:
system (
Atoms) – The ASE Atoms object representing the systemdyn (dynamics object) – The initialized MD dynamics object
steps (
int) – The number of MD steps to runpace (
int) – The interval (in steps) for logging and saving checkpointslog_filename (
str) – The filename for the MD log filetrajfile (
str) – The filename for the trajectory filedir_name (
str) – The directory where log and trajectory files will be savedname (
str) – A label for the simulation (e.g., thermostat type)temp_start (float | None)
temp_end (float | None)
- sparc.src.ase_md.ExecuteMlpDynamics(system, dyn, steps, pace, log_filename, trajfile, dir_name, distance_metrics, name, epot_threshold, temp_start=None, temp_end=None, restart=False)[source]
Run a machine learning potential MD simulation with safety checks.
Supports restart from checkpoint saved inside dir_name.
- Parameters:
system (
Atoms) – The ASE Atoms object representing the systemdyn (dynamics object) – The initialized MD dynamics object
steps (
int) – The number of MD steps to runpace (
int) – The interval (in steps) for logginglog_filename (
str) – The filename for the MD log filetrajfile (
str) – The filename for the trajectory filedir_name (
str) – The directory where log, trajectory and checkpoint files are saveddistance_metrics (
Optional[List[Dict]]) – Metrics used to check physical limits during simulationname (
str) – A label for the simulation (e.g., thermostat type)epot_threshold (
float) – Maximum allowed energy deviation from reference (eV)restart (
bool) – If True, resume from checkpoint in dir_name (default: False)temp_start (float | None)
temp_end (float | None)
- sparc.src.ase_md.LangevinNVT(atoms, timestep=1, temperature=300, friction=0.01, restart=False)[source]
Set up a Langevin thermostat for NVT MD simulation.
- Parameters:
atoms (
Atoms) – The ASE Atoms object representing the systemtimestep (
float) – The simulation timestep in femtoseconds (default: 1 fs)temperature (
float) – The target temperature in Kelvin (default: 300 K)friction (
float) – The friction coefficient in fs^-1 (default: 0.01)restart (
bool) – If True, restart from checkpoint (default: False)
- Returns:
The initialized Langevin dynamics object
- Return type:
dynamics
- sparc.src.ase_md.NPT(system, timestep, temperature, tau_t, pressure, tau_p, compressibility=None, restart=False, **kwargs)[source]
Set up NPT dynamics using ASE’s NPTBerendsen integrator.
- Parameters:
system (
Atoms) – The ASE Atoms object representing the systemtimestep (
float) – MD integration timestep in femtosecondstemperature (
float) – Target temperature in Kelvintau_t (
float) – Thermostat time constant in femtosecondspressure (
float) – Target pressure in bartau_p (
float) – Barostat time constant in femtosecondscompressibility (
Optional[float]) – Isothermal compressibility in 1/bar. If None, uses the default for Cu (~7.1e-7 bar⁻¹)restart (
Optional[bool]) – If True, restart from checkpoint (default: False)**kwargs (dict) – Extra ASE MD options
- Returns:
The initialized NPTBerendsen dynamics object
- Return type:
dynamics
- sparc.src.ase_md.NVE(system, timestep, restart=False)[source]
Set up NVE dynamics using ASE’s VelocityVerlet integrator.
- sparc.src.ase_md.NoseNVT(atoms, timestep=1, temperature=300, tdamp=10, restart=False)[source]
Set up a Nose-Hoover chain NVT thermostat for MD simulation.
- Parameters:
atoms (
Atoms) – The ASE Atoms object representing the systemtimestep (
float) – The simulation timestep in femtoseconds (default: 1 fs)temperature (
float) – The target temperature in Kelvin (default: 300 K)tdamp (
float) – The damping time for the thermostat in femtoseconds (default: 10 fs)restart (
bool) – If True, restart from checkpoint (default: False)
- Returns:
The initialized NoseHooverChainNVT dynamics object
- Return type:
dynamics
- sparc.src.ase_md.TemperatureRamp(dyn, atoms, current_step, total_steps, temp_start=None, temp_end=None, ensemble='NVT')[source]
Apply linear temperature ramping to MD dynamics object.
- Implements VASP-style temperature ramping:
T(t) = T_start + (T_end - T_start) * (t / t_total)
Temperature ramping is only applicable to NVT ensemble.
- Parameters:
dyn (dynamics object) – MD dynamics object (NoseHooverChainNVT, Langevin, NPTBerendsen)
atoms (
Atoms) – ASE Atoms objectcurrent_step (
int) – Current MD steptotal_steps (
int) – Total MD stepstemp_start (
Optional[float]) – Starting temperature (K). If None, no ramping appliedtemp_end (
Optional[float]) – Ending temperature (K). If None, no ramping appliedensemble (
str) – MD ensemble (‘NVT’, ‘NPT’, ‘NVE’). Default: ‘NVT’
- sparc.src.ase_md.initialize_dynamics(atoms, dyn_class, timestep, restart=False, **kwargs)[source]
Initialize MD dynamics, optionally restarting from a checkpoint.
This is a unified initialization function that handles velocity initialization, checkpoint loading, and dynamics object creation for all ensemble types.
- Parameters:
atoms (
Atoms) – The ASE Atoms object representing the systemdyn_class (type) – The MD dynamics class (NoseHooverChainNVT, Langevin, NPTBerendsen, VelocityVerlet)
timestep (
float) – The simulation timestep in ASE internal units (already converted)restart (
bool) – If True, restart from checkpoint (default: False)**kwargs (dict) – Additional parameters for specific ensemble types
- Returns:
The initialized dynamics object
- Return type:
dyn
References
For more information on ASE, visit: https://wiki.fysik.dtu.dk/ase/