Molecular Dynamics

class quippy.dynamicalsystem.DynamicalSystem(name=dynamicalsystem)[source]

Bases: DynamicalSystem

Defined at DynamicalSystem.fpp lines 195-228

Attributes:
atoms

Element atoms ftype=type(atoms) pytype=Atoms

avg_temp

Element avg_temp ftype=real(dp) pytype=float

avg_time

Element avg_time ftype=real(dp) pytype=float

cur_temp

Element cur_temp ftype=real(dp) pytype=float

dt

Element dt ftype=real(dp) pytype=float

dw

Element dw ftype=real(dp) pytype=float

ekin

Element ekin ftype=real(dp) pytype=float

epot

Element epot ftype=real(dp) pytype=float

ext_energy

Element ext_energy ftype=real(dp) pytype=float

group_lookup

Element group_lookup ftype=integer pytype=int

initialised

Element initialised ftype=logical pytype=bool

n

Element n ftype=integer pytype=int

nconstraints

Element nconstraints ftype=integer pytype=int

ndof

Element ndof ftype=integer pytype=int

nrestraints

Element nrestraints ftype=integer pytype=int

nrigid

Element nrigid ftype=integer pytype=int

nsteps

Element nsteps ftype=integer pytype=int

print_thermostat_temps

Element print_thermostat_temps ftype=logical pytype=bool

random_seed

Element random_seed ftype=integer pytype=int

t

Element t ftype=real(dp) pytype=float

thermostat_dw

Element thermostat_dw ftype=real(dp) pytype=float

thermostat_work

Element thermostat_work ftype=real(dp) pytype=float

wkin

Element wkin ftype=real(dp) pytype=float

work

Element work ftype=real(dp) pytype=float

Methods

add_atoms(*args, **kwargs)

Defined at DynamicalSystem.fpp lines 235-236

advance_verlet(self, dt, f[, virial, e, ...)

Defined at DynamicalSystem.fpp lines 1867-1887

advance_verlet1(self, dt[, virial, parallel, ...)

Defined at DynamicalSystem.fpp lines 1547-1714

advance_verlet2(self, dt, f[, virial, e, ...)

Defined at DynamicalSystem.fpp lines 1740-1857

constrain_atom_plane(self, i, plane_n[, d, ...)

Defined at DynamicalSystem.fpp lines 2225-2248

constrain_bondanglecos(self, i, j, k[, c, ...)

Defined at DynamicalSystem.fpp lines 2005-2038

constrain_bondlength(self, i, j, d[, di, t0, ...)

Defined at DynamicalSystem.fpp lines 2041-2082

constrain_bondlength_dev_pow(self, i, j, p, ...)

Defined at DynamicalSystem.fpp lines 2114-2155

constrain_bondlength_diff(self, i, j, k, d[, ...)

Defined at DynamicalSystem.fpp lines 2161-2195

constrain_bondlength_sq(self, i, j[, d, ...)

Defined at DynamicalSystem.fpp lines 2085-2111

constrain_gap_energy(self, d[, restraint_k, ...)

Defined at DynamicalSystem.fpp lines 2201-2222

constrain_struct_factor_like_i(self, z, q, ...)

Defined at DynamicalSystem.fpp lines 2347-2392

constrain_struct_factor_like_mag(self, z, q, ...)

Defined at DynamicalSystem.fpp lines 2251-2296

constrain_struct_factor_like_r(self, z, q, ...)

Defined at DynamicalSystem.fpp lines 2299-2344

disable_damping(self)

Defined at DynamicalSystem.fpp lines 994-996

ds_add_constraint(self, atoms, func, data[, ...)

Defined at DynamicalSystem.fpp lines 2396-2468

ds_add_thermostat(self, type_bn, t[, gamma, ...)

Defined at DynamicalSystem.fpp lines 842-889

ds_add_thermostats(self, type_bn, n[, t, ...)

Defined at DynamicalSystem.fpp lines 896-942

ds_amend_constraint(self, constraint, func, ...)

Defined at DynamicalSystem.fpp lines 2474-2479

ds_angular_momentum(self[, origin, indices])

Defined at DynamicalSystem.fpp lines 1041-1046

ds_kinetic_energy(self[, local_ke, error])

Defined at DynamicalSystem.fpp lines 1135-1143

ds_kinetic_virial(self[, error])

Defined at DynamicalSystem.fpp lines 1194-1201

ds_momentum(self[, indices])

Defined at DynamicalSystem.fpp lines 1005-1009

ds_print(self[, file])

Defined at DynamicalSystem.fpp lines 1970-1994

ds_print_status(self[, label, epot, ...)

Defined at DynamicalSystem.fpp lines 1899-1963

ds_print_thermostats(self)

Defined at DynamicalSystem.fpp lines 1965-1967

ds_remove_thermostat(self, index_bn)

Defined at DynamicalSystem.fpp lines 837-840

ds_restore_state(self, from_)

Defined at DynamicalSystem.fpp lines 604-625

ds_save_state(self, from_[, error])

Defined at DynamicalSystem.fpp lines 564-596

ds_set_barostat(self, type_bn, p_ext, ...)

Defined at DynamicalSystem.fpp lines 814-830

ds_update_thermostat(self[, t, p, i])

Defined at DynamicalSystem.fpp lines 949-957

enable_damping(self, damp_time)

Defined at DynamicalSystem.fpp lines 988-992

get_damping_time(self)

Defined at DynamicalSystem.fpp lines 979-984

is_damping_enabled(self)

Defined at DynamicalSystem.fpp lines 970-977

n_thermostat(self)

Defined at DynamicalSystem.fpp lines 832-835

remove_atoms(*args, **kwargs)

Defined at DynamicalSystem.fpp lines 240-241

rescale_velo(self, temp[, mass_weighted, zero_l])

Defined at DynamicalSystem.fpp lines 1342-1370

temperature(self[, property, value, ...)

Defined at DynamicalSystem.fpp lines 1270-1324

thermostat_temperatures(self, temps)

Defined at DynamicalSystem.fpp lines 1246-1258

zero_momentum(self[, indices])

Defined at DynamicalSystem.fpp lines 1422-1434

add_atoms(*args, **kwargs)[source]

Defined at DynamicalSystem.fpp lines 235-236

Overloaded interface containing the following procedures:

_ds_add_atom_single _ds_add_atom_multiple

Add one or more atoms to this DynamicalSystem. Equivalent to ‘Atoms%add_atoms’, but also appends the number of degrees of freedom correctly.

advance_verlet(self, dt, f[, virial, e, parallel, store_constraint_force, do_calc_dists, error])[source]

Defined at DynamicalSystem.fpp lines 1867-1887

Parameters:
dsDynamicalsystem
dtfloat
ffloat array
virialfloat array
efloat
parallelbool
store_constraint_forcebool
do_calc_distsbool
errorint
Calls’ advance_verlet2’ followed by ‘advance_verlet1’. Outside this routine the
velocities will be half-stepped. This allows a simpler MD loop:
> for n in range(n_steps):
> pot.calc(ds.atoms, force=True, energy=True)
> ds.advance_verlet(dt, ds.atoms.force)
> ds.print_status(epot=ds.atoms.energy)
> if n % connect_interval == 0:
> ds.atoms.calc_connect()
advance_verlet1(self, dt[, virial, parallel, store_constraint_force, do_calc_dists, error])[source]

Defined at DynamicalSystem.fpp lines 1547-1714

Parameters:
thisDynamicalsystem
dtfloat
virialfloat array
parallelbool
store_constraint_forcebool
do_calc_distsbool
errorint
Advance the velocities by half the time-step ‘dt’ and the
positions by a full time-step. A typical MD loop should
resemble the following code(Python example, Fortran is similar):
> ds.atoms.calc_connect()
> for n in range(n_steps):
> ds.advance_verlet1(dt)
> pot.calc(ds.atoms, force=True, energy=True)
> ds.advance_verlet2(dt, ds.atoms.force)
> ds.print_status(epot=ds.atoms.energy)
> if n % connect_interval == 0:
> ds.atoms.calc_connect()
advance_verlet2(self, dt, f[, virial, e, parallel, store_constraint_force, error])[source]

Defined at DynamicalSystem.fpp lines 1740-1857

Parameters:
thisDynamicalsystem
dtfloat
ffloat array
virialfloat array
efloat
parallelbool
store_constraint_forcebool
errorint
Advances the velocities by the second half time-step
constrain_atom_plane(self, i, plane_n[, d, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2225-2248

Parameters:
thisDynamicalsystem
iint
plane_nfloat array
dfloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain an atom to lie in a particluar plane
constrain_bondanglecos(self, i, j, k[, c, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2005-2038

Parameters:
thisDynamicalsystem
iint
jint
kint
cfloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain the bond angle cosine between atoms i, j, and k
constrain_bondlength(self, i, j, d[, di, t0, tau, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2041-2082

Parameters:
thisDynamicalsystem
iint
jint
dfloat
difloat
t0float
taufloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain the bond between atoms i and j
constrain_bondlength_dev_pow(self, i, j, p, d[, di, t0, tau, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2114-2155

Parameters:
thisDynamicalsystem
iint
jint
pfloat
dfloat
difloat
t0float
taufloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain the bond between atoms i and j
constrain_bondlength_diff(self, i, j, k, d[, di, t0, tau, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2161-2195

Parameters:
thisDynamicalsystem
iint
jint
kint
dfloat
difloat
t0float
taufloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain the difference of bond length between atoms i–j and j–k
constrain_bondlength_sq(self, i, j[, d, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2085-2111

Parameters:
thisDynamicalsystem
iint
jint
dfloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain the bond between atoms i and j
constrain_gap_energy(self, d[, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2201-2222

Parameters:
thisDynamicalsystem
dfloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain the energy gap of two resonance structures
constrain_struct_factor_like_i(self, z, q, sf[, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2347-2392

Parameters:
thisDynamicalsystem
zint
qfloat array
sffloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain an atom to lie in a particluar plane
constrain_struct_factor_like_mag(self, z, q, sf[, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2251-2296

Parameters:
thisDynamicalsystem
zint
qfloat array
sffloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain an atom to lie in a particluar plane
constrain_struct_factor_like_r(self, z, q, sf[, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2299-2344

Parameters:
thisDynamicalsystem
zint
qfloat array
sffloat
restraint_kfloat
boundint
tolfloat
print_summarybool
Constrain an atom to lie in a particluar plane
disable_damping(self)[source]

Defined at DynamicalSystem.fpp lines 994-996

Parameters:
thisDynamicalsystem
ds_add_constraint(self, atoms, func, data[, update_ndof, restraint_k, bound, tol, print_summary])[source]

Defined at DynamicalSystem.fpp lines 2396-2468

Parameters:
thisDynamicalsystem
atomsint array
funcint
datafloat array
update_ndofbool
restraint_kfloat
boundint
tolfloat
print_summarybool
Add a constraint to the DynamicalSystem and reduce the number of degrees of freedom,
unless ‘update_Ndof’ is present and false.
ds_add_thermostat(self, type_bn, t[, gamma, q, tau, tau_cell, p, bulk_modulus_estimate, cell_oscillation_time, nhl_tau, nhl_mu, massive, region_i])[source]

Defined at DynamicalSystem.fpp lines 842-889

Parameters:
thisDynamicalsystem
type_bnint
tfloat
gammafloat
qfloat
taufloat
tau_cellfloat
pfloat
bulk_modulus_estimatefloat
cell_oscillation_timefloat
nhl_taufloat
nhl_mufloat
massivebool
region_iint
Add a new thermostat to this DynamicalSystem. ‘type’ should
be one of the following thermostat types:
* ‘THERMOSTAT_NONE’
* ‘THERMOSTAT_LANGEVIN’
* ‘THERMOSTAT_NOSE_HOOVER’
* ‘THERMOSTAT_NOSE_HOOVER_LANGEVIN’
* ‘THERMOSTAT_LANGEVIN_NPT’
* ‘THERMOSTAT_LANGEVIN_PR’
* ‘THERMOSTAT_NPH_ANDERSEN’
* ‘THERMOSTAT_NPH_PR’
* ‘THERMOSTAT_LANGEVIN_OU’
* ‘THERMOSTAT_LANGEVIN_NPT_NB’
‘T’ is the target temperature. ‘Q’ is the Nose-Hoover coupling constant. Only one
of ‘tau’ or ‘gamma’ should be given. ‘p’ is the external
pressure for the case of Langevin NPT.
ds_add_thermostats(self, type_bn, n[, t, t_a, gamma, gamma_a, q, q_a, tau, tau_a, region_i])[source]

Defined at DynamicalSystem.fpp lines 896-942

Parameters:
thisDynamicalsystem
type_bnint
nint
tfloat
t_afloat array
gammafloat
gamma_afloat array
qfloat
q_afloat array
taufloat
tau_afloat array
region_iint
ds_amend_constraint(self, constraint, func, data[, k])[source]

Defined at DynamicalSystem.fpp lines 2474-2479

Parameters:
thisDynamicalsystem
constraintint
funcint
datafloat array
kfloat
Replace a constraint involving some atoms with a different constraint involving
the same atoms
ds_angular_momentum(self[, origin, indices])[source]

Defined at DynamicalSystem.fpp lines 1041-1046

Parameters:
thisDynamicalsystem
originfloat array
indicesint array
Returns:
lfloat array
Return the angular momentum of all the atoms in this DynamicalSystem, defined by
$mathbf{L} = sum_{i} mathbf{r_i} imes mathbf{v_i}$.
ds_kinetic_energy(self[, local_ke, error])[source]

Defined at DynamicalSystem.fpp lines 1135-1143

Parameters:
thisDynamicalsystem
local_kebool
errorint
Returns:
kefloat

Return the total kinetic energy $E_k = sum_{i}

rac{1}{2} m v^2$
ds_kinetic_virial(self[, error])[source]

Defined at DynamicalSystem.fpp lines 1194-1201

Parameters:
thisDynamicalsystem
errorint
Returns:
kvfloat array

Return the total kinetic virial $w_ij = sum_{k}

rac{1}{2} m v_i v_j$
ds_momentum(self[, indices])[source]

Defined at DynamicalSystem.fpp lines 1005-1009

Parameters:
thisDynamicalsystem
indicesint array
Returns:
pfloat array
Return the total momentum $mathbf{p} = sum_i mathbf{m_i} mathbf{v_i}$.
Optionally only include the contribution of a subset of atoms.
ds_print(self[, file])[source]

Defined at DynamicalSystem.fpp lines 1970-1994

Parameters:
thisDynamicalsystem
fileInoutput
Print lots of information about this DynamicalSystem in text format.
ds_print_status(self[, label, epot, instantaneous, file, error])[source]

Defined at DynamicalSystem.fpp lines 1899-1963

Parameters:
thisDynamicalsystem
labelstr
epotfloat
instantaneousbool
fileInoutput
errorint
Print a status line showing the current time, temperature, the mean temperature
the total energy and the total momentum for this DynamicalSystem. If present, the optional
‘label’ parameter should be a one character label for the log lines and is printed
in the first column of the output. ‘epot’ should be the potential energy
if this is available. ‘instantaneous’ has the same meaning as in ‘temperature’ routine.
ds_print_thermostats(self)[source]

Defined at DynamicalSystem.fpp lines 1965-1967

Parameters:
thisDynamicalsystem
ds_remove_thermostat(self, index_bn)[source]

Defined at DynamicalSystem.fpp lines 837-840

Parameters:
thisDynamicalsystem
index_bnint
ds_restore_state(self, from_)[source]

Defined at DynamicalSystem.fpp lines 604-625

Parameters:
toDynamicalsystem
from_Dynamicalsystem
Restore a DynamicalSystem to a previously saved state.
Only scalar members and ‘ds%atoms’ (minus ‘ds%atoms%connect’)
are copied back; ‘to’ should be a properly initialised
DynamicalSystem object. The saved state of the random
number generator is also restored. ‘calc_dists()’ is
called on the restored atoms object.
ds_save_state(self, from_[, error])[source]

Defined at DynamicalSystem.fpp lines 564-596

Parameters:
toDynamicalsystem
from_Dynamicalsystem
errorint
Save the state of a DynamicalSystem. The output object
cannot be used as an initialised DynamicalSystem since
connectivity and group information is not copied to save
memory. Only scalar members and the ‘ds%atoms’ object
(minus ‘ds%atoms%connect’) are copied. The current
state of the random number generator is also saved.
ds_set_barostat(self, type_bn, p_ext, hydrostatic_strain, diagonal_strain, finite_strain_formulation, tau_epsilon[, w_epsilon, t, w_epsilon_factor, thermalise])[source]

Defined at DynamicalSystem.fpp lines 814-830

Parameters:
thisDynamicalsystem
type_bnint
p_extfloat
hydrostatic_strainbool
diagonal_strainbool
finite_strain_formulationbool
tau_epsilonfloat
w_epsilonfloat
tfloat
w_epsilon_factorfloat
thermalisebool
ds_update_thermostat(self[, t, p, i])[source]

Defined at DynamicalSystem.fpp lines 949-957

Parameters:
thisDynamicalsystem
tfloat
pfloat
iint
enable_damping(self, damp_time)[source]

Defined at DynamicalSystem.fpp lines 988-992

Parameters:
thisDynamicalsystem
damp_timefloat
Enable damping, with damping time set to `damp_time`. Only atoms
flagged in the `damp_mask` property will be affected.
get_damping_time(self)[source]

Defined at DynamicalSystem.fpp lines 979-984

Parameters:
thisDynamicalsystem
Returns:
get_damping_timefloat
is_damping_enabled(self)[source]

Defined at DynamicalSystem.fpp lines 970-977

Parameters:
thisDynamicalsystem
Returns:
is_damping_enabledbool
n_thermostat(self)[source]

Defined at DynamicalSystem.fpp lines 832-835

Parameters:
thisDynamicalsystem
Returns:
n_thermostatint
remove_atoms(*args, **kwargs)[source]

Defined at DynamicalSystem.fpp lines 240-241

Overloaded interface containing the following procedures:

_ds_remove_atom_single _ds_remove_atom_multiple

Remove one or more atoms from this DynamicalSystem. Equivalent of ‘Atoms%remove_atoms’, but also amends the number of degrees of freedom correctly.

rescale_velo(self, temp[, mass_weighted, zero_l])[source]

Defined at DynamicalSystem.fpp lines 1342-1370

Parameters:
thisDynamicalsystem
tempfloat
mass_weightedbool
zero_lbool
Rescale the atomic velocities to temperature ‘temp’. If the
current temperature is zero, we first randomise the velocites.
If ‘mass_weighted’ is true, then the velocites are weighted by
$1/sqrt{m}$. Linear momentum is zeroed automatically. If
‘zero_l’ is true then the angular momentum is also zeroed.
temperature(self[, property, value, include_all, instantaneous, error])[source]

Defined at DynamicalSystem.fpp lines 1270-1324

Parameters:
thisDynamicalsystem
propertystr
valueint
include_allbool
instantaneousbool
errorint
Returns:
temperaturefloat

Return the temperature, assuming each degree of freedom contributes $

rac{1}{2}kT$. By default only moving and thermostatted atoms are

included — this can be overriden by setting ‘include_all’ to true. ‘region’ can be used to restrict the calculation to a particular thermostat region. ‘instantaneous’ controls whether the calculation should be carried out using the current values of the velocities and masses, or whether to return the value at the last Verlet step (the latter is the default).

thermostat_temperatures(self, temps)[source]

Defined at DynamicalSystem.fpp lines 1246-1258

Parameters:
thisDynamicalsystem
tempsfloat array
zero_momentum(self[, indices])[source]

Defined at DynamicalSystem.fpp lines 1422-1434

Parameters:
thisDynamicalsystem
indicesint array
Change velocities to those that the system would have in the zero momentum frame.
Optionalally zero the total momentum of a subset of atoms, specified by ‘indices’.
property atoms

Element atoms ftype=type(atoms) pytype=Atoms

Defined at DynamicalSystem.fpp line 221

property avg_temp

Element avg_temp ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 205

Time-averaged temperature

property avg_time

Element avg_time ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 207

Averaging time, in fs

property cur_temp

Element cur_temp ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 206

Current temperature

property dt

Element dt ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 204

Last time step

property dw

Element dw ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 208

Increment of work done this time step

property ekin

Element ekin ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 211

Current kinetic energy

property epot

Element epot ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 210

Total potential energy

property ext_energy

Element ext_energy ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 213

Extended energy

property group_lookup

Element group_lookup ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 219

Stores which group atom $i$ is in

property initialised

Element initialised ftype=logical pytype=bool

Defined at DynamicalSystem.fpp line 216

property n

Element n ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 197

Number of atoms

property nconstraints

Element nconstraints ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 200

Number of constraints

property ndof

Element ndof ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 202

Number of degrees of freedom

property nrestraints

Element nrestraints ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 201

Number of restraints

property nrigid

Element nrigid ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 199

Number of rigid bodies

property nsteps

Element nsteps ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 198

Number of integration steps

property print_thermostat_temps

Element print_thermostat_temps ftype=logical pytype=bool

Defined at DynamicalSystem.fpp line 228

property random_seed

Element random_seed ftype=integer pytype=int

Defined at DynamicalSystem.fpp line 217

RNG seed, used by ‘ds_save_state’ and ‘ds_restore_state’ only.

property t

Element t ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 203

Time

property thermostat_dw

Element thermostat_dw ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 214

Increment of work done by thermostat

property thermostat_work

Element thermostat_work ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 215

Total work done by thermostat

property wkin

Element wkin ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 212

Current kinetic contribution to the virial

property work

Element work ftype=real(dp) pytype=float

Defined at DynamicalSystem.fpp line 209

Total work done

class quippy.dynamicalsystem.Dynamics(atoms, timestep, trajectory, trajectoryinterval=10, initialtemperature=None, logfile='-', loginterval=1, loglabel='D')[source]

Bases: Dynamics

Wrapper around DynamicalSystem integrator compatible with ASE ase.md.MolecularDynamics interface.

Note

This class uses ASE units for time (and hence velocities and momenta) which differ from the femtoseconds used in QUIP. When using this class, all times should be given in ASE time units.

Parameters:
  • atoms – quippy or ASE Atoms instance. Note that if atoms is not a quippy Atoms object, a copy will be made.

  • timestep – in ASE time units (multiply by ase.units.fs to convert from femtoseconds)

  • trajectory – output file to which to write the trajectory. Can be a string to create a new output object.

  • trajectoryinterval – interval at which to write frames

  • initialtemperature – if not None, rescale initial velocities to a specified temperature, given in K.

  • logfile – filename or open InOutput object to write log lines to. Use special filename "-" for stdout (default).

  • loginterval – interval at which to write log lines

Attributes:
average_temperature

Average temperature

averaging_time

Averaging time used for average temperature and kinetic energy

damping

Get or set the damping time constant in fs.

number_of_constraints

Get number of constraints

number_of_degrees_of_freedom

Get number of degrees of freedom in system, including QUIP constraints

number_of_restraints

Get number of restraints

number_of_rigid_bodies

Get number of rigid_bodies

state

Save or restore current state of this dynamical system

temperature

Get or set the current temperature

time

Time in ASE units (Note: NOT the same as femtoseconds)

timestep

Set timestep, in ASE time units

Methods

add_thermostat(type, T[, gamma, Q, tau, ...])

Add a new thermostat to this Dynamics, and return its index.

attach(function[, interval])

Attach callback function.

converged()

MD is 'converged' when number of maximum steps is reached.

get_number_of_thermostats()

Return the number of active thermostats

get_thermostat_temperatures()

Return the current temperatures of all thermostated regions

insert_observer(function[, position, interval])

Insert an observer.

irun()

Run dynamics algorithm as generator.

log(*args)

a dummy function as placeholder for a real logger, e.g.

print_thermostats()

Print a table of the current thermostats in this system

remove_thermostat(index)

Remove thermostat number index.

run([steps])

Run dynamics forwards for steps steps.

set_barostat(type, p_ext, ...[, W_epsilon, ...])

Add a barostat to this dynamical system

set_temperature(temperature)

Randomise velocities to a target temperature (given in K)

step(forces)

Advance dynamics by one time-step.

update_barostat(p, T)

Update target pressure p or temperature T for NPT barostat

update_thermostat([T, p, index])

Change the target temperature T or presure p for thermostat index.

add_thermostat(type, T, gamma=None, Q=None, tau=None, tau_cell=None, p=None, NHL_tau=None, NHL_mu=None, massive=None)[source]

Add a new thermostat to this Dynamics, and return its index.

By default, all thermostats apply to the whole system. The ‘thermostat_region’ property can be used to control which thermostats affect which atoms. It should be set to the index of the required thermostat for each atom.

Parameters:
  • type – string or integer giving thermostat type. The following types are supported:: THERMOSTAT_NONE, THERMOSTAT_LANGEVIN, THERMOSTAT_NOSE_HOOVER, THERMOSTAT_NOSE_HOOVER_LANGEVIN, THERMOSTAT_LANGEVIN_NPT, THERMOSTAT_LANGEVIN_PR, THERMOSTAT_NPH_ANDERSEN, THERMOSTAT_NPH_PR, THERMOSTAT_LANGEVIN_OU, THERMOSTAT_LANGEVIN_NPT_NB, THERMOSTAT_ALL_PURPOSE.

  • T – target temperature for this thermostat, in K.

  • gamma – decay constant, in units of 1/fs

  • tau – time constant, in units of fs. tau == 1/gamma, and only one of the two should be specified.

  • tau_cell – time constant for the cell degrees of freedom (for variable cell thermostats only).

  • p – target pressure (for variable cell thermostats)

  • NHL_tau – time constant for Nose-Hoover-Langevein thermostats

  • NHL_mu – thermostat mass Nose-Hoover-Langevin thermostats

  • massive – set to True to enable massive Nose-Hoover-Langevin

attach(function, interval=1, *args, **kwargs)[source]

Attach callback function.

At every interval steps, call function with arguments args and keyword arguments kwargs.

converged()[source]

MD is ‘converged’ when number of maximum steps is reached.

get_number_of_thermostats()[source]

Return the number of active thermostats

Excludes thermostat with index zero, which is used for damping.

get_thermostat_temperatures()[source]

Return the current temperatures of all thermostated regions

insert_observer(function, position=0, interval=1, *args, **kwargs)[source]

Insert an observer.

irun()[source]

Run dynamics algorithm as generator. This allows, e.g., to easily run two optimizers or MD thermostats at the same time.

Examples: >>> opt1 = BFGS(atoms) >>> opt2 = BFGS(StrainFilter(atoms)).irun() >>> for _ in opt2: >>> opt1.run()

log(*args)[source]

a dummy function as placeholder for a real logger, e.g. in Optimizer

print_thermostats()[source]

Print a table of the current thermostats in this system

remove_thermostat(index)[source]

Remove thermostat number index.

index should be in range 1 <= index <= dyn.get_number_of_thermostats()

run(steps=50)[source]

Run dynamics forwards for steps steps.

set_barostat(type, p_ext, hydrostatic_strain, diagonal_strain, finite_strain_formulation, tau_epsilon, W_epsilon=None, T=None, W_epsilon_factor=None)[source]

Add a barostat to this dynamical system

Langevin NPT barostat based on fluctuation/dissipation, as originally implemented for hydrostatic strain using EOM from Quigley, D. and Probert, M.I.J., J. Chem. Phys. 120 11432, subsequently rediscretized by Noam Bernstein based on ideas from and discussion with B. Leimkuhler.

Modified for arbitrary strain based on (but not exactly using) formulation in E. Tadmor and R. Miller Modeling Materials: Continuum, Atomistic and Multiscale Techniques (Cambridge University Press, 2011). Chap 9.5. Their definition of stress (factors of F, F^T, J, and F^-T) for finite strain is optionally used, but certain terms in EOM are excluded, and their discretization is entirely ignored.

Parameters:
  • type – The type of barostat to be used. Currently only BAROSTAT_HOOVER_LANGEVIN is supported.

  • p_ext – External pressure in QUIP units eV/A^3

  • hystrostatic_strain – Set to True to constrain a hydrostatic strain tensor

  • diagonal_strain – Set to True to constrain a diagonal strain tensor

  • finite_strain_formulation – If True, use finite instead of infinitessimal strain formulation

  • tau_epsilon – time constant for Langevin part (in fs)

  • W_epsilon – fictious cell mass used in Langevin NPT

  • T – target temperature for Langevin part

  • W_epsilon_factor – Scaling factor for fictious cell mass

set_temperature(temperature)[source]

Randomise velocities to a target temperature (given in K)

step(forces)[source]

Advance dynamics by one time-step.

Returns forces at the end of the time-step, suitable for passing to next step()

update_barostat(p, T)[source]

Update target pressure p or temperature T for NPT barostat

update_thermostat(T=None, p=None, index=1)[source]

Change the target temperature T or presure p for thermostat index.

property average_temperature

Average temperature

property averaging_time

Averaging time used for average temperature and kinetic energy

property damping

Get or set the damping time constant in fs. Set to None to disable damping. By default damping applies to all atoms, but can be selectively enabled with the ‘damp_mask’ property.

property number_of_constraints

Get number of constraints

property number_of_degrees_of_freedom

Get number of degrees of freedom in system, including QUIP constraints

property number_of_restraints

Get number of restraints

property number_of_rigid_bodies

Get number of rigid_bodies

property state

Save or restore current state of this dynamical system

property temperature

Get or set the current temperature

property time

Time in ASE units (Note: NOT the same as femtoseconds)

property timestep

Set timestep, in ASE time units