AMS Molecular Dynamics PLAMS jobs

The following special Jobs exist to simplify running MD simulations (and continuing/restarting MD simulations):

  • AMSMDJob. A general class for which you can set up most common MD options using arguments to the constructor.

  • AMSNVEJob for NVE simulations

  • AMSNVTJob for NVT simulations

  • AMSNPTJob for NPT simulations

  • AMSMDScanDensityJob for running MD deformation simulations while isotropically scaling the density (from AMS2025 a similar input option is available in the AMS Driver, MolecularDynamics%Deformation%TargetDensity)

  • AMSNVESpawnerJob is a special MultiJob that runs several NVE simulations with initial velocities taken from evenly spaced frames in a previous job.

Some default values are different from AMS. For example, the checkpoint frequency is set to a higher number, and the thermostat constant tau is automatically set to 400 times the timestep, by default.

Always check the input by calling job.get_input() before running a job.

In general, the jobs and results classes inherit from AMSJob and AMSResults.

The NVE, NVT, and NPT classes simply remove the below options if they are set:

allowed input block

AMSMDJob

NVE

NVT

NPT

Thermostat

yes

no

yes

yes

Barostat

yes

no

no

yes

Nanoreactor

yes

no

no

no

Deformation

yes

no

no

no

The following jobs help with the postanalysis:

  • AMSMSDJob for calculating mean square displacement (MSD)

  • AMSRDFJob for calculating radial distribution functions (RDF)

  • AMSVACFJob for calculating velocity autocorrelation functions (VACF)

AMSMDJob API

class AMSMDJob(velocities=None, timestep=None, samplingfreq=None, nsteps=None, checkpointfrequency=None, engineresultsfreq=None, writevelocities=None, writebonds=None, writemolecules=None, writecharges=None, writeenginegradients=None, calcpressure=None, molecule=None, temperature=None, thermostat=None, tau=None, thermostat_region=None, pressure=None, barostat=None, barostat_tau=None, scale=None, equal=None, constantvolume=None, binlog_time=None, binlog_pressuretensor=None, binlog_dipolemoment=None, _enforce_thermostat=False, _enforce_barostat=False, **kwargs)[source]
molecule: Molecule

The initial structure

name: str

The name of the job

settings: Settings

Settings for the job. You should normally not populate neither settings.input.ams.MolecularDynamics nor settings.input.ams.Task

velocities: float or AMSJob or str (path/to/ams.rkf) or 2-tuple (path/to/ams.rkf, frame-number)

If float, it is taken as the temperature. If AMSJob or str, the velocities are taken from the EndVelocities section of the corresponding ams.rkf file. If 2-tuple, the first item must be a path to an ams.rkf, and the second item an integer specifying the frame number - the velocities are then read from History%Velocities(framenumber).

timestep: float

Timestep

samplingfreq: int

Sampling frequency

nsteps: int

Length of simulation

Trajectory options:

checkpointfrequency: int

How frequently to write MDStep*.rkf checkpoint files

writevelocitiesbool

Whether to save velocities to ams.rkf (needed for example to restart from individual frames or to calculate velocity autocorrelation functions)

writebonds: bool

Whether to save bonds to ams.rkf

writemolecules: bool

Whether to write molecules to ams.rkf

writeenginegradients: bool

Whether to save engine gradients (negative of atomic forces) to ams.rkf

Thermostat (NVT, NPT) options:

thermostat: str

‘Berendsen’ or ‘NHC’

tau: float

Thermostat time constant (fs)

temperature: float or tuple of floats

Temperature (K). If a tuple/list of floats, the Thermostat.Duration option will be set to evenly divided intervals.

thermostat_region: str

Region for which to apply the thermostat

Barostat (NPT) options:

barostat: str

‘Berendsen’ or ‘MTK’

barostat_tau: float

Barostat time constant (fs)

pressure: float

Barostat pressure (pascal)

equal: str

‘XYZ’ etc.

scale: str

‘XYZ’ etc.

constantvolume: bool

Constant volume?

Other options:

calcpressure: bool

Whether to calculate pressure for each frame.

binlog_time: bool

Whether to log the time at every timestep in the BinLog section on ams.rkf

binlog_pressuretensor: bool

Whether to log the pressure tensor at every timestep in the BinLog section on ams.rkf

get_velocities_from(other_job, frame=None, update_molecule=True)[source]

Function to update the InitialVelocities block in self. It is normally not needed, instead use the e.g. AMSNVEJob.restart_from() function.

This function can be called in prerun() methods for MultiJobs

classmethod restart_from(other_job, frame=None, settings=None, use_prerun=False, **kwargs)[source]
other_job: AMSJob

The job to restart from.

frame: int

Which frame to read the structure and velocities from. If None, the final structure and end velocities will be used (section Molecule and MDResults%EndVelocities)

settings: Settings

Settings that override any other settings. All settings from other_job (e.g. the engine settings) are inherited by default but they can be overridden here.

use_prerun: bool

If True, the molecule and velocities will only be read from other_job inside the prerun() method. Set this to True to prevent PLAMS from waiting for other_job to finish as soon as the new job is defined.

kwargs: many options

See the docstring for AMSMDJob.

AMSNVEJob API

This class uses the same arguments as AMSMDJob.

class AMSNVEJob(**kwargs)[source]

A class for running NVE MD simulations

AMSNVTJob API

This class uses the same arguments as AMSMDJob.

class AMSNVTJob(**kwargs)[source]

A class for running NVT MD simulations

AMSNPTJob API

This class uses the same arguments as AMSMDJob.

AMSNPTResults inherits from AMSResults.

class AMSNPTJob(**kwargs)[source]

A class for running NPT MD simulations

class AMSNPTResults(*args, **kwargs)[source]
get_equilibrated_molecule(equilibration_fraction=0.667, return_index=False)[source]

Discards the first equilibration_fraction of the trajectory. Calculates the average density of the rest. Returns the molecule with the closest density to the average density among the remaining trajectory.

AMSNVESpawnerJob API

class AMSNVESpawnerJob(previous_job, n_nve=1, name='nvespawnerjob', **kwargs)[source]

A class for running multiple NVE simulations with initial structures/velocities taken from an NVT trajectory. The NVT trajectory must contain the velocities!

__init__(previous_job, n_nve=1, name='nvespawnerjob', **kwargs)[source]
previous_job: AMSJob

An AMSJob with an MD trajectory. Must contain velocities (WriteVelocities Yes). Note that the trajectory should have been equilibrated before it starts.

n_nveint

The number of NVE simulations to spawn

All other settings can be set as for an AMSNVEJob (e.g. nsteps).

AMSMDScanDensityJob API

class AMSMDScanDensityJob(molecule, scan_density_upper=1.4, startstep=None, **kwargs)[source]

A class for scanning the density using MD Deformations

AMSRDFJob API

class AMSRDFJob(previous_job, atom_indices=None, atom_indices_to=None, rmin=0.5, rmax=6.0, rstep=0.1, **kwargs)[source]
__init__(previous_job, atom_indices=None, atom_indices_to=None, rmin=0.5, rmax=6.0, rstep=0.1, **kwargs)[source]
previous_job: AMSJob

AMSJob with finished MD trajectory.

atom_indices: List[int]

Atom indices (starting with 1). If None, calculate RDF from all atoms.

atom_indices_to: List[int]

Atom indices (starting with 1). If None, calculate RDF to all atoms.

rmin: float

Minimum distance (angstrom)

rmax: float

Maximum distance (angstrom)

rstep: float

Bin width for the histogram (angstrom)

class AMSRDFResults(job)[source]

Results class for AMSRDFJob

get_rdf()[source]

Returns a 2-tuple r, rdf.

r: numpy array of float (angstrom) rdf: numpy array of float

AMSMSDJob API

class AMSMSDJob(previous_job, max_correlation_time_fs=10000, start_time_fit_fs=2000, atom_indices=None, **kwargs)[source]

A convenient class wrapping around the trajectory analysis MSD tool

__init__(previous_job, max_correlation_time_fs=10000, start_time_fit_fs=2000, atom_indices=None, **kwargs)[source]
previous_job: AMSJob

An AMSJob with an MD trajectory. Note that the trajectory should have been equilibrated before it starts.

max_correlation_time_fs: float

Maximum correlation time in femtoseconds

start_time_fit_fsfloat

Smallest correlation time for the linear fit

atom_indicesList[int]

If None, use all atoms. Otherwise use the specified atom indices (starting with 1)

kwargs: dict

Other options to AMSAnalysisJob

class AMSMSDResults(job)[source]

Results class for AMSMSDJob

get_msd()[source]

returns time [fs], msd [ang^2]

get_linear_fit(start_time_fit_fs=None, end_time_fit_fs=None)[source]

Fits the MSD between start_time_fit_fs and end_time_fit_fs

Returns a 3-tuple LinRegress.result, fit_x (fs), fit_y (ang^2)

result.slope is given in ang^2/fs

get_diffusion_coefficient(start_time_fit_fs=None, end_time_fit_fs=None)[source]

Returns D in m^2/s

AMSVACFJob API

class AMSVACFJob(previous_job, max_correlation_time_fs=5000, max_freq=5000, atom_indices=None, **kwargs)[source]

A class for calculating the velocity autocorrelation function and its power spectrum

__init__(previous_job, max_correlation_time_fs=5000, max_freq=5000, atom_indices=None, **kwargs)[source]
previous_job: AMSJob

An AMSJob with an MD trajectory. Note that the trajectory should have been equilibrated before it starts.

max_correlation_time_fs: float

Maximum correlation time in femtoseconds

max_freq: float

The maximum frequency for the power spectrum in cm^-1

atom_indices: List[int]

Atom indices (starting with 1). If None, use all atoms.

class AMSVACFResults(job)[source]

Results class for AMSVACFJob

get_vacf()[source]

Extract the velocity autocorrelation function vs. time.

Returns a 2-tuple time, y with time in fs.

get_power_spectrum(max_freq=None)[source]

Calculate the power spectrum as the Fourier transform of the velocity autocorrelation function.

max_freq: float, optional

The maximum frequency in cm^-1.

Returns: A 2-tuple freq, y with freq in cm^-1 and y unitless.

AMSViscosityFromBinLogJob API

class AMSViscosityFromBinLogJob(previous_job, **kwargs)[source]

A convenient class wrapping around the trajectory analysis ViscosityFromBinLog tool. Only runs with default input options (max correlation time 10% of the total trajectory).

__init__(previous_job, **kwargs)[source]
previous_job: AMSJob

An AMSJob with an MD trajectory. Note that the trajectory should have been equilibrated before it starts. It must have been run with the BinLog%PressureTensor option set.

kwargs: dict

Other options to AMSAnalysisJob

class AMSViscosityFromBinLogResults(job)[source]

Results class for AMSViscosityFromBinLogJob

get_viscosity_integral()[source]

Extract the running viscosity integral.

Returns a 2-tuple with 1D numpy arrays time, viscosity_integral, with time in fs and viscosity_integral in Pa*s.

get_double_exponential_fit()[source]

Perform a double exponential fit to the viscosity integral.

The fitted function is of the form

A * (lam * (1 - np.exp(-x / tau1)) + (1 - lam) * (1 - np.exp(-x / tau2)))

where A is the limiting value in the infinite time limit.

Returns: a 3-tuple popt, time, prediction, where

  • popt is a 4-tuple containing A, lam, tau1, and tau2,

  • time is in fs, and

  • prediction is the value of the above function vs time.