MDMC.MD package

Subpackages

Submodules

MDMC.MD.container module

Module for the AtomContainer class

class MDMC.MD.container.AtomContainer[source]

Bases: ABC

A collection of Atom objects

The AtomContainer can be indexed and returns an Atom. However indexing cannot be used for setting or deleting an Atom.

atoms

A list of the Atom objects that belong to the AtomContainer

Type:

list

abstract property atoms: list[Atom]

returns: A list of the Atom objects belonging to the AtomContainer :rtype: list

MDMC.MD.interaction_functions module

A module for storing atomic interaction functions

Contains class InteractionFunction from which all interaction function classes must derive. All functions describing atomic interactions must be added to this module in order to be called by a universe. If needed, the interaction function classes can be extended to contain actual function definitions.

Contains class Parameter, which defines the name and value of each parameter which belongs to an InteractionFunction, and whether the parameter is fixed, has constraints or is tied.

Contains filters for filtering list of parameters based on a predicate.

class MDMC.MD.interaction_functions.Buckingham(A: float, B: float, C: float)[source]

Bases: InteractionFunction

The Buckingham potential (in units of kJ mol^-1) for the interaction of 2 atoms at distance r (in Ang) has the form:

\[{\Phi _{12}(r)=A\exp \left(-Br\right)-{\frac {C}{r^{6}}}}\]

The first and second terms represent a repulsion and attraction respectively, and so all parameters \(A\), \(B\) and \(C\) should be positive in order to be physically valid.

Parameters:
  • A (float) – The Buckingham parameter A in units of kJ mol^-1

  • B (float) – The Buckingham parameter B in units of Ang^-1

  • C (float) – The Buckingham parameter C in units of Ang^6 kJ mol^-1

Examples

The following creates a Dispersion interaction with a Buckingham functional form:

buck = Buckingham(3.65e-18, 6.71, 6.94e-22)
O_disp = Dispersion(universe, (O.atom_type, O.atom_type), function=buck)
class MDMC.MD.interaction_functions.Coulomb(charge: float)[source]

Bases: InteractionFunction

Coulomb interaction for charged particles:

\[E = \frac{Cq_{i}q_{j}}{r}\]
Parameters:

charge (float) – The charge in units of e

Examples

The following creates a Coulombic interaction with a Coulomb functional form:

O = Atom('O')
coul = Coulomb(-0.8476)
O_coulombic = Coulombic(atoms=O, cutoff=10., function=coul)

As Coulomb is the default functional form of Coulombic interactions when an Atom is created, this is equivalent to setting the charge on an Atom:

O = Atom('O', charge=-0.8476)
class MDMC.MD.interaction_functions.HarmonicPotential(equilibrium_state: float, potential_strength: float, **settings: dict)[source]

Bases: InteractionFunction

Harmonic potential for bond stretching, and angular and improper dihedral vibration, with the form:

\[E = K(r-r_0)^2\]

As HarmonicPotential can be used with several different Interaction types, the Interaction type must be specified, so that the correct units can be assigned to the equilibrium_state and potential_strength parameters.

Parameters:
  • equilibrium_state (float) – The equilibrium state of the object in either Ang or degrees, depending on the interaction_type passed.

  • potential_strength (float) – The potential strength in units of kJ mol^-1 Ang^-2 (linear) or kJ mol^-1 rad^-2 (angular/improper), depending on the interaction_type passed.

  • **settings

    interaction_typestr

    A str specifying either 'bond', 'angle' or 'improper'. This assigns the correct units to equilibrium_state and potential_strength. This keyword must be passed.

Raises:
  • ValueError – The interaction_type must be 'bond', 'angle', or 'improper'

  • TypeError – An interaction_type of 'bond', 'angle', or 'improper' must be passed

Examples

The following creates a HarmonicPotential for a Bond interaction:

hp = HarmonicPotential(1., 4637., interaction_type='bond')

The following creates a HarmonicPotential for a BondAngle interaction:

hp = HarmonicPotential(109.47, 383., interaction_type='angle')

The following creates a HarmonicPotential for a DihedralAngle interaction with improper==True (i.e. an improper dihedral):

hp = HarmonicPotential(180., 20.92, interaction_type='improper')
class MDMC.MD.interaction_functions.InteractionFunction(val_dict: dict)[source]

Bases: object

Base class for interaction functions, which can be user supplied

Parameters:

val_dict (dict) – name:value pairs. Currently this must be ordered alphabetically. value must either be an object with a value and a unit (e.g. a UnitFloat object), or a (float, str) tuple, where float is the value and str is the unit.

property name: str

Get the name of the class of the InteractionFunction

Returns:

The class name

Return type:

str

property parameters: Parameters

Get or set the interaction function’s parameters

Returns:

A Parameters object containing each Parameter

Return type:

Parameters

property parameters_values: ndarray

Get the values for all Parameters objects

Returns:

A NumPy array of values for all Parameter

Return type:

numpy.ndarray

set_parameters_interactions(interaction: Interaction) None[source]

Sets the parent Interaction for all Parameters objects

Parameters:

interaction (Interaction) – An Interaction to set as the parent of all the Parameters

class MDMC.MD.interaction_functions.LennardJones(epsilon: float, sigma: float, **settings: dict)[source]

Bases: InteractionFunction

Dispersive Lennard-Jones interaction with the form:

\[E = 4{\epsilon}[(\frac{\sigma}{r})^{12} - (\frac{\sigma}{r})^6)] \qquad r < r_c\]
Parameters:
  • epsilon (float) – The LJ epsilon value in units of kJ mol^-1

  • sigma (float) – The LJ sigma value in units of Ang

  • **settings

    cutofffloat

    The distance in Ang at which the potential is cutoff

    long_range_solverstr

    The long range solver, either 'PPPM', 'PME', or 'E' for Particle-Particle Particle-Mesh, Particle Mesh Ewald, or Ewald solvers

Examples

The following creates a Dispersion interaction with a LennardJones functional form:

lj = LennardJones(0.6502, 3.166)
O_disp = Disperion(universe, (O.atom_type, O.atom_type), function=lj)
class MDMC.MD.interaction_functions.Periodic(K1: float, n1: float, d1: float, *parameters: float)[source]

Bases: InteractionFunction

Periodic potential for proper and improper dihedral interactions, with the form:

\[E = \sum_{i=1,m}K_i[1.0+\cos(n_i\phi-d_i)]\]

where phi is angle between the ijk and jkl planes (where i, j, k, and l are the four atoms of the dihedral).

Parameters:
  • K1 (float) – The K parameter (energy) of the first order term, in units of kJ mol^-1

  • n1 (int) – The n parameter of the first order term, which is unitless. This must be a non-negative int.

  • d1 (float) – The d parameter (angle) of the first order term, in units of deg

  • *parameters – K, n, and d parameters for higher order terms. These must be ordered K2, n2, d2, K3, n3, d3, K4, n4, d4 etc. The types and units of these parameters are the same as the corresponding first order parameters listed above.

Examples

The following creates a first order Periodic for a DihedralAngle interaction with improper==True (i.e. an improper dihedral):

periodic = Periodic(87.864, 2, 180.)
improper = DihedralAngle(atoms=[C, H1, H2, O], improper=True,
                         function=periodic)

The following creates a third order Periodic for a DihedralAngle interaction with improper==False (i.e. a proper dihedral), with K1=3.53548, K2=-4.02501 and K3=2.98319:

periodic = Periodic(3.53548, 1, 0., -4.02501, 2, 180., 2.98319, 3, 0.)
proper = DihedralAngle(atoms=[N, C1, C2, C3], improper=False,
                       function=periodic)
MDMC.MD.interaction_functions.inter_func_decorator(*parameter_units) Callable[source]

Decorates a method to add units to all positional and any relevant keyword arguments

Designed for adding units to parameters of __init__ method for subclasses of InteractionFunction.

Parameters:

*parameter_units (tuple) – One or more tuple where the first element is a str giving the keyword name of an expected parameter. The second element is a str or Unit, where each str (or Unit) is a unit which is applied to the corresponding value passed to the decorated method. If one of the values is unitless, pass None at the corresponding index in parameter_units.

Examples

The following adds units of 'Ang' to parameter alpha, units of 's' to the parameter beta, and units of 'atm' to the parameter gamma:

@inter_func_decorator(('alpha', 'Ang'), ('beta', 's'), ('gamma', 'Pa'))
def __init__(self, alpha, beta, gamma):
    ...

If one of the parameters is unitless, this can be set with None (in which case the returned type will be the same as the original value i.e. a UnitFloat or UnitNDArray will not be created). So to set epsilon as unitless:

@inter_func_decorator(('delta', 'arb'), ('epsilon', None), ('gamma', 'deg'))
def __init__(self, delta, epsilon, gamma):
    ...

MDMC.MD.interactions module

Module in which all interactions between structural units are defined.

Interaction is the abstract base class from which all interactions have to be derived..

class MDMC.MD.interactions.Bond(*atom_tuples: tuple, **settings: dict)[source]

Bases: ConstrainableMixin, BondedInteraction

A bond between any two atoms. Requires exactly two atoms in each atom_tuple.

Parameters:
  • atom_tuples (list) – A list of tuple. Each tuple contains Atom which are bonded together.

  • **settings

is_equivalent(other) bool[source]

Checks for equivalence between two BondedInteraction``s, specifically if they apply to the same ``atom_types, have the same function describing the interaction, and have the same constrained setting.

Parameters:

other (BondedInteraction) – The object to compare against.

Return type:

bool

class MDMC.MD.interactions.BondAngle(*atom_tuples: tuple, **settings: dict)[source]

Bases: ConstrainableMixin, BondedInteraction

A bond angle between any two bonds

Requires three Atom objects (rotation around central atom) in each atom_tuple. The atoms are ordered i, j, k, where j is the central atom. So:

BondAngle(i, j, k) == BondAngle(k, j, i)
Parameters:
  • atom_tuples (list) – A list of tuple. Each tuple contains Atom which are bonded together. For three or more Atom objects, the order of the Atom objects within each tuple is important.

  • **settings

is_equivalent(other: BondedInteraction) bool[source]

Checks for equivalence between two BondedInteraction``s, specifically if they apply to the same ``atom_types, have the same function describing the interaction, and have the same constrained setting.

Parameters:

other (BondedInteraction) – The object to compare against.

Return type:

bool

class MDMC.MD.interactions.BondedInteraction(*atom_tuples: list[tuple], **settings: dict)[source]

Bases: Interaction

Base class for bonded interactions

Parameters:
  • atom_tuples (list) – A list of tuple. Each tuple contains Atom objects which are bonded together. For three or more Atom objects, the order of the Atom objects within each tuple is important.

  • **settings

    n_atoms (int)

    The number of atoms to which this BondedInteraction applies, for example 2 for a Bond.

Examples

For a single bonded interactions which applies to H1, O1, and H2:

BondedInteraction(H1, O1, H2)

For two bonded interactions of the same BondedInteraction type, one applied to H1, O1 and H2 Atom objects, and the other applied to H3, O2 and H4 Atom objects:

BondedInteraction((H1, O1, H2), (H3, O2, H4))

Whereas the above examples are both specifying a H-O-H ordered BondedInteraction, the following specifies a H-H-O BondedInteraction:

BondAngle(H1, H2, O)
add_atoms(*atoms: Atom, **settings: dict) None[source]

Add atoms which are all involved in one example of this interaction

Parameters:
  • *atoms – one or more Atom objects

  • **settings

    from_structure (bool)

    If add_atoms has been called from a Structure

Raises:

ValueError – If this BondedInteraction has already been applied to one or more of the atoms

property atom_types: Set[int | None]

Get the set of all atom_type``s that this ``BondedInteraction corresponds to, including None if appropriate.

Returns:

A set of all (unique) ``atom_type``s.

Return type:

Set[Union[int, NoneType]]

property atoms: list[tuple[Atom]]

Get or set the atoms on which the Coulombic interaction is applied

Returns:

A list of tuple containing one or more Atom. Each tuple contains all of the atoms involved in one example of the interaction. For example a BondAngle interaction each tuple would contain 3 or 4 atoms.

Return type:

list

Raises:

TypeError – If a list of tuple is not set

element_list() list | None[source]

Get a list of the elements for which the BondedInteraction applies

Returns:

The elements for which the BondedInteraction applies

Return type:

list

is_equivalent(other) bool[source]

Checks for equivalence between two BondedInteraction``s, specifically if they apply to the same ``atom_types and the same function describing the interaction.

Parameters:

other (BondedInteraction) – The object to compare against.

Return type:

bool

property universe: Universe | None

Get the Universe to which the BondedInteraction belongs

Returns:

The Universe to which the BondedInteraction belongs or None

Return type:

Universe

class MDMC.MD.interactions.ConstrainableMixin(*atom_tuples: tuple, **settings: dict)[source]

Bases: object

A mixin class enabling classes inheriting from BondedInteraction to be constrained

These constraints are then applied by a constraint algorithm (e.g. SHAKE), which is specified in the Universe to which the BondedInteraction belongs.

Parameters:
  • atom_tuples (list) – A list of tuple. Each tuple contains Atom objects which are bonded together. For three or more Atom objects, the order of the Atom objects within each tuple is important.

  • **settings

constrained

Specifying whether the object is constrained

Type:

bool

class MDMC.MD.interactions.Coulombic(universe: Universe = None, **settings: dict)[source]

Bases: NonBondedInteraction

A non-bonded coulombic interaction - either normal or modified Coulomb

Parameters:
  • universe (Universe, optional) – The Universe in which the Coulombic exists. Default is None. Must be passed as a parameter if atom_types if passed.

  • **settings

    charge (float)

    The charge parameter of the Coulombic interaction, in units of e. If this argument is passed, the interaction_function of this Coulombic is set to Coulomb with this float as its Parameter. Passing charge will overwrite any other interaction_functions that are set, i.e. it makes function parameter redundant

    atoms (list)

    Atom objects to which the Coulombic applies. If specifying the atoms, universe does not need to be passed as a parameter.

    atom_typeslist of int

    int for each atom_type for which the NonBondedInteraction applies. If specifying the atom_types, the universe must be passed as a parameter and the atoms for which the atom_types are specified must exist in Universe. See the example above in the ‘charge’ section.

Raises:
  • TypeError – If one or more atom_types are passed but no universe is passed

  • TypeError – If neither atom_types or atoms have been passed

  • TypeError – If both atom_types and atoms have been passed

Examples

Upon initializing an Atom object and adding it to a Universe:

O = Atom('O', atom_type=1)
universe = Universe(10.0)
universe.add_structure('O')

The following initializations of Coulombic are equivalent:

O_coulombic = Coulombic(universe, atom_types=[O.atom_type],
                        charge=-0.84)
O_coulombic = Coulombic(universe, atom_types=[O.atom_type],
                        function=Coulomb(-0.84))

If atoms is passed then a Universe does not need to be passed:

O_coulombic = Coulombic(atoms=[O], charge=-0.84)
property atom_types: list

Get the atom types for which the Coulombic applies

Returns:

All atom types to which the Coulombic applies. If the interaction was initialized with atoms, all atom_types of the atoms to which the Coulombic was applied are returned; HOWEVER THE COULOMBIC INTERACTION IS NOT APPLIED TO ALL ATOMS OF THESE atom_types, ONLY THE ATOMS IN self.atoms

Return type:

list

property atoms: list[Atom]

Get the atoms on which the Coulombic interaction is applied

Returns:

A list of Atom on which the Coulombic is applied

Return type:

list

element_list() list[source]

Get a list of the elements for which the Coulombic interaction applies

Returns:

The elements for which the Coulombic interaction applies

Return type:

list

class MDMC.MD.interactions.DihedralAngle(*atom_tuples: tuple, **settings: dict)[source]

Bases: BondedInteraction

A dihedral angle between any two sets of three atoms, ijk and jkl.

Dihedral angles can be both proper and improper, where the angle between the two planes of ijk and jkl is fixed for improper dihedrals.

The atoms of a proper DihedralAngle are ordered i, j, k, l, where j and k are the two central atoms. So:

DihedralAngle(i, j, k, l) == DihedralAngle(l, k, j, i)

The atoms of an improper DihedralAngle are ordered i, j, k, l, where i is the central atom to which j, k, and l are all connected. So:

(DihedralAngle(i, j, k, l, improper=True)
== DihedralAngle(i, j, l, k, improper=True)
== DihedralAngle(i, l, k, j, improper=True)
== DihedralAngle(i, l, j, k, improper=True))
== DihedralAngle(i, k, j, l, improper=True))
== DihedralAngle(i, k, l, j, improper=True))
Parameters:
  • atom_tuples (list) – A list of tuple. Each tuple contains four Atom objects which are bonded together by the DihedralAngle, in the order specified.

  • **settings

    improper (bool)

    Whether the DihedralAngle is improper or not.

improper

Whether the DihedralAngle is improper or not, which affects the InteractionFunction which can be set for this DihedralAngle. By default this is set to False i.e. the interaction is a proper dihedral.

Type:

bool

class MDMC.MD.interactions.Dispersion(universe: Universe, *atom_types: int, **settings: dict)[source]

Bases: NonBondedInteraction

A non-bonded dispersive interaction - either LJ or Buckingham

Parameters:
  • universe (Universe) – The Universe in which the NonBondedInteraction exists

  • *atom_typesint for each atom type for which the NonBondedInteraction applies

  • **settings

    cutoff (float)

    The distance in Ang at which the interaction potential is truncated

    vdw_tail_correction (bool)

    Specifies if the tail correction to the energy and pressure should be applied. This only affects the simulation dynamics if the simulation is being performed with constant pressure.

Raises:
  • TypeErroratom_types must be iterable

  • ValueErrorDispersion should only be specified as existing between pairs of atom_types

  • TypeError – Each atom_type must be int

property atom_types

Get the atom types for which the NonBondedInteraction applies

Returns:

A list of int for the atom_types

Return type:

list

property atoms: list[tuple[list[Atom]]]

Get the atoms on which the Dispersion is applied

Returns:

A list of two tuple, where each tuple contains a list of Atom. Every Atom in the first tuple has a dispersion interaction with every Atom in the second tuple (excluding self interactions). This is the complete list of possible dispersion interactions, i.e. it is only exactly correct if no cutoff has been specified.

Return type:

list

element_list() list[source]

Get a list of the elements for which the Interaction applies

Returns:

The elements for which the Interaction applies

Return type:

list

is_equivalent(other) bool[source]

Checks for equivalence between two Dispersion``s, specifically if they apply to the same ``atom_types, have the same cuttoff, the same function describing the interaction and vdw_tail_correction setting.

Parameters:

other (Dispersion) – The object to compare against.

Return type:

bool

class MDMC.MD.interactions.Interaction(**settings: dict)[source]

Bases: ABC

Base class for interactions, both bonded, non-bonded and constraints

Each different type of interaction should have an Interaction object. This object contains a list of the Atom (or `Atom pairs, triplets or quadruplets, depending on the type of interaction) for which this Interaction applies. For example, an oxygen Coulombic interaction would contain a list of tuple where each tuple contains a different O Atom, and a hydrogen-oxygen Bond interaction would contain a list of tuple where each tuple contains a different H and O pair. Interaction objects can be sliced to return a sublist of the tuple.

When an Atom is passed to an Interaction, the Interaction is also added to the Atom.

Parameters:

**settings

function (InteractionFunction)

A class of interaction function (e.g. HarmonicPotential)

abstract property atoms: Atom | list[Atom]

Get the atoms on which the Interaction is applied

abstract element_list() list[source]

Get a list of the elements for which the Interaction applies

Returns:

The elements for which the Interaction applies

Return type:

list

element_tuple() tuple[source]

A tuple of elements for which the Interaction applies

Returns:

The elements for which the Interaction applies

Return type:

tuple

property function: InteractionFunction

Get or set the InteractionFunction of the Interaction

Returns:

The interaction function of the Interaction

Return type:

InteractionFunction

property function_name: str | None

Get the name of the InteractionFunction belonging to the Interaction

Returns:

The name of the InteractionFunction, or None if no InteractionFunction has been set

Return type:

str

property parameters: Parameters

Get the Parameter objects belonging to the InteractionFunction belonging to the Interaction

Returns:

A Parameters object containing each Parameter

Return type:

Parameters

sorted_element_list() list[source]

Sort the list of elements for which the Interaction applies

Returns:

The elements for which the Interaction applies, sorted alphabetically

Return type:

list

abstract property universe: Universe

Get the Universe to which the Interaction belongs

Returns:

The Universe to which the Interaction belongs or None

Return type:

Universe

class MDMC.MD.interactions.NonBondedInteraction(universe, *atom_types, **settings)[source]

Bases: Interaction

Base class for non-bonded interactions

Parameters:
  • universe (Universe) – The Universe in which the NonBondedInteraction exists

  • *atom_typesint for each atom_type for which the NonBondedInteraction applies

  • **settings

    cutoff (float)

    The distance in Ang at which the interaction potential is truncated

abstract property atom_types: list[int]

Get the atom types for which the NonBondedInteraction applies

Returns:

A list of int for the atom_types

Return type:

list

property cutoff: float

Get or set the distance in Ang at which the interaction potential is truncated

Returns:

The distance in Ang of the cutoff

Return type:

float

is_equivalent(other) bool[source]

Checks for equivalence between two NonBondedInteraction``s, specifically if they apply to the same ``atom_types, have the same cuttoff and the same function describing the interaction.

Parameters:

other (NonBondedInteraction) – The object to compare against.

Return type:

bool

property universe: Universe

Get or set the Universe to which the NonBondedInteraction belongs

Returns:

The Universe to which the NonBondedInteraction belongs or None

Return type:

Universe

MDMC.MD.parameters module

A module for the Parameter and Parameters classes

Parameter defines the name and value of each force field parameter, and whether it is fixed, has constraints or is tied.

Parameters inherits from lists and implements a number of methods for filterting a sequence of Parameter objects.

class MDMC.MD.parameters.Parameter(value, name, fixed=False, constraints=None, **settings)[source]

Bases: object

A force field parameter which can be fixed or constrained within limits

The value of a parameter cannot be set if fixed==True.

Parameters:
  • value (float) – The value of the parameter.

  • name (str) – The name of the parameter.

  • fixed (bool) – Whether or not the value can be changed.

  • constraints (tuple) – The closed range of the Parameter.value, (lower, upper). constraints must have the same units as value.

  • **settings

    unit (str)

    The unit. If this is not provided then the unit will be taken from the object passed as value.

property constraints
property interactions: list

Get or append to the parent Interaction objects for this Parameter

Returns:

All parent Interaction objects

Return type:

list

Raises:
  • ValueError – If an added interaction name is not consistent with existing interaction names

  • ValueError – If an added Interaction has a function name not consistent with the function names of an existing Interaction

set_tie(parameter: Parameter, expr: str) None[source]

This ties the Parameter.value to the value of another Parameter

Parameters:
  • parameter (Parameter) – The Parameter to tie to

  • expr (str) – A mathematical expression

Examples

To set the Parameter.value to p1.value * 2:

>>> Parameter.set_tie(p1, "* 2")
property tie: float | None

Get the value of a the Parameter that this Parameter is tied to

Returns:

The value of the tied Parameter

Return type:

float

property tied: bool

Get whether this Parameter is tied

Returns:

True if this Parameter is tied to another Parameter, else False

Return type:

bool

static validate_value(value: float, constraints: tuple) None[source]

Validates the Parameter.value by testing if it is within the constraints

Parameters:
  • values (float) – The value of the Parameter

  • constraints (tuple) – A 2-tuple of the lower and upper constraints respectively.

Raises:

ValueError – If the value is not within the constraints

property value: float

Get or set the value of the Parameter

The value will not be changed if it is fixed or tied, or if it is set outside the bounds of constraints

Returns:

The value of the Parameter, including if the Parameter is tied

Return type:

float

Warns:
  • warnings.warn – If the Parameter is fixed.

  • warnings.warn – If the Parameter is tied.

class MDMC.MD.parameters.Parameters(init_parameters: Parameter | list[Parameter] | None = None)[source]

Bases: dict

A dict-like object where every element is a Parameter indexed by name, which contains a number of helper methods for filtering.

Although Parameters is a dict, it should be treated like a list when writing to it; i.e. initialise it using a list and use append to add to it. These parameters can then be accessed by their name as a key.

In short; Parameters writes like a list and reads like a dict.

Parameters:

init_parameters (Parameter or list of ``Parameter``s, optional, default None) – The initial Parameter objects that the Parameters object contains.

array

An alphabetically-sorted numpy array of the ``Parameter``s stored in this object.

Type:

np.ndarray

append(parameters: list[Parameter] | Parameter) None[source]

Appends a Parameter or list of ``Parameter``s to the dict, with the parameter name as its key.

Parameters:

parameters (Parameter or list of ``Parameter``s) – The parameter(s) to be added to the dict.

property as_array: ndarray

The parameters in the object as a sorted numpy array.

Returns:

An alphabetically-sorted array of parameter values in the object.

Return type:

np.ndarray

filter(predicate: Callable[[Parameter], bool]) Parameters[source]

Filters using a predicate

Parameters:

predicate (function) – A function that returns a bool which takes a Parameter as an argument.

Returns:

The Parameter objects which meet the condition of the predicate

Return type:

Parameters

filter_atom_attribute(attribute: str, value: str | float) Parameters[source]

Filters based on the attribute of Atom objects which have each Parameter applied to them

Parameters:
  • attribute (str) – An attribute of an Atom. Attributes to match to must be either float or str.

  • value (str, float) – The value of the Atom attribute.

Returns:

The Parameter objects which are applied to an Atom object which has the specified value of the specified attribute

Return type:

Parameters

filter_function(function_name: str) Parameters[source]

Filters based on the name of the InteractionFunction of each Parameter

Parameters:

function_name (str) – The name of the InteractionFunction of Parameter objects to return, for example 'LennardJones' or 'HarmonicPotential'.

Returns:

The Parameter objects which have a function with the specified function_name

Return type:

Parameters

filter_interaction(interaction_name: str) Parameters[source]

Filters based on the name of the Interaction of each Parameter

Parameters:

interaction_name (str) – The name of the Interaction of Parameter objects to return, for example 'Bond'.

Returns:

The Parameter objects which have an Interaction with the specified interaction_name

Return type:

Parameters

filter_name(name: str) Parameters[source]

Filters by name

Parameters:

name (str) – The name of the Parameter objects to return.

Returns:

The Parameter objects with name

Return type:

Parameters

filter_structure(structure_name: str) Parameters[source]

Filters based on the name of the Structure to which each Parameter applies

Parameters:

structure_name (str) – The name of a Structure.

Returns:

The Parameter objects which are applied to a Structure which has the specified structure_name

Return type:

Parameters

filter_value(comparison: str, value: float) Parameters[source]

Filters by value

Parameters:
  • comparison (str) – A str representing a comparison operator, '>', '<', '>=', '<=', '==', '!='.

  • value (float) – A float with which Parameter values are compared, using the comparison operator.

Returns:

The Parameter objects which return a True when their values are compared with value using the comparison operator

Return type:

Parameters

log_parameters() None[source]

Logs all Parameters by ID

MDMC.MD.simulation module

Module for setting up and running the simulation

Classes for the simulation box, minimizer and integrator.

class MDMC.MD.simulation.ConstraintAlgorithm(accuracy: float, max_iterations: int)[source]

Bases: object

Class describing the algorithm and parameters which are applied to constrain BondedInteraction objects

Parameters:
  • accuracy (float) – The accuracy (tolerance) of the applied constraints

  • max_iterations (int) – The maximum number of iterations that can be used when calculating the additional force that is required to constrain the atoms to satisfy the constraints on the bonded interactions

accuracy

The accuracy (tolerance) of the applied constraints

Type:

float

property max_iterations: int

Get or set the maximum number of iterations that can be used when calculating the additional force that is required to constrain the atoms to satisfy the constraints on the bonded interactions

Returns:

The maximum number of iterations

Return type:

int

property name: str

Get the name of the class

Returns:

The name of the class

Return type:

str

class MDMC.MD.simulation.Ewald(**settings: dict)[source]

Bases: KSpaceSolver

Holds the parameters that are required for the Ewald solver to be applied to both/either the electrostatic and/or dispersion interactions

Parameters:

**settings

accuracy (float)

The relative RMS error in per-atom forces

class MDMC.MD.simulation.KSpaceSolver(**settings: dict)[source]

Bases: object

Class describing the k-space solver that is applied to electrostatic and/or dispersion interactions

Different MDEngine require different parameters to be specified for a k-space solver to be used. These parameters are specified in settings.

Parameters:

**settings

accuracy (float)

The relative RMS error in per-atom forces

accuracy

The relative RMS error in per-atom forces

Type:

float

property name

Get the name of the class

Returns:

The name of the class

Return type:

str

class MDMC.MD.simulation.PPPM(**settings: dict)[source]

Bases: KSpaceSolver

Holds the parameters that are required for the PPPM solver to be applied to both/either the electrostatic and/or dispersion interactions

Parameters:

**settings

accuracy (float)

The relative RMS error in per-atom forces

class MDMC.MD.simulation.Rattle(accuracy: float, max_iterations: int)[source]

Bases: ConstraintAlgorithm

Holds the parameters which are required for the RATTLE algorithm to be applied to the constrained interactions

Parameters:
  • accuracy (float) – The accuracy (tolerance) of the applied constraints

  • max_iterations (int) – The maximum number of iterations that can be used when calculating the additional force that is required to constrain the atoms to satisfy the constraints on the bonded interactions

class MDMC.MD.simulation.Shake(accuracy: float, max_iterations: int)[source]

Bases: ConstraintAlgorithm

Holds the parameters which are required for the SHAKE algorithm to be applied to the constrained interactions

Parameters:
  • accuracy (float) – The accuracy (tolerance) of the applied constraints

  • max_iterations (int) – The maximum number of iterations that can be used when calculating the additional force that is required to constrain the atoms to satisfy the constraints on the bonded interactions

class MDMC.MD.simulation.Simulation(universe: Universe, traj_step: int, time_step: float = 1.0, engine: str = 'lammps', **settings)[source]

Bases: object

Sets up the molecular dynamics engine and parameters for how it should run

An ensemble is defined by whether a thermostat or barostat are present

Parameters:
  • universe (Universe) – The Universe on which the simulation is performed.

  • traj_step (int) – How many steps the simulation should take between dumping each CompactTrajectory frame. Along with time_step determines the time separation of calculated variables such as energy.

  • time_step (float, optional) – Simulation timestep in fs. Default is 1.

  • engine (str, optional) – The MDEngine used for the simulation. Default is 'lammps'.

  • **settings

    temperature (float)

    Simulation temperature in K. Note that velocities of atoms in the MD engine are set based on the temperature. If all atoms in the universe have 0 velocity, then velocities for the MD engine will be randomly chosen from a uniform distribution and then scaled to give the correct temperature. If one or more velocities are set, then these will be scaled to give the correct temperature. In either case, only the velocities on the engine, not the universe, are affected.

    integrator (str)

    Simulation time integrator.

    lj_options (dict)

    option:value pairs for Lennard-Jones interactions.

    es_options (dict)

    option:value pairs for electrostatic interactions.

    thermostat (str)

    The name of the thermostat e.g. Nose-Hoover.

    barostat (str)

    The name of the barostat e.g. Nose-Hoover.

    pressure (float)

    Simulation pressure in Pa. This is required if a barostat is passed.

    verbose (str)

    If the output of the class instantiation should be reported, default to True.

universe

The Universe on which the simulation is performed.

Type:

Universe

traj_step

How many steps the simulation should take between dumping each CompactTrajectory frame. Along with time_step determines the time separation of calculated variables such as energy.

Type:

int

engine

A subclass of MDEngine which provides the interface to the MD library. Default is 'lammps'.

Type:

MDEngine, optional

settings

The settings passed to the Simulation. See the Parameters section for details.

Type:

dict

time_step
trajectory
auto_equilibrate(variables: list[str] = ['temp', 'pe'], eq_step: int = 10, window_size: int = 100, tolerance: float = 0.01) Tuple[int, dict][source]

Equilibrate until the specified list of variables have stabilised. Uses the KPSS stationarity test to determine whether the variables are stationary in a given window.

Parameters:
  • variables (list[str], default ['temp', 'pe']) – The MD engine variables we use to monitor stability.

  • eq_step (int, default 10) – The number of equilibration steps between each stability check.

  • window_size (int, default 100) – The size of the rolling window of stability checks. The simulation is considered stable if it has been stationary for eq_step*window_size equilibration steps.

  • tolerance (float, 0.01) – The p-value used for the KPSS test.

Returns:

  • int – The number of equilibration steps performed.

  • vals_dict – The value of each variable per step, for graphing if desired.

minimize(n_steps: int, minimize_every: int = 10, verbose: bool = False, output_log: str = None, work_dir: str = None, **settings: dict) None[source]

Performs an MD run intertwined with periodic structure relaxation. This way after a local minimum is found, the system is taken out of the minimum to explore a larger volume of the parameter space.

Parameters:
  • n_steps (int) – Total number of the MD run steps

  • minimize_every (int, optional) – Number of MD steps between two consecutive minimizations

  • verbose (bool, optional) – Whether to print statements when the minimization has been started and completed (including the number of minimization steps and time taken). Default is False.

  • output_log (str, optional) – Log file for the MD engine to write to. Default is None.

  • work_dir (str, optional) – Working directory for the MD engine to write to. Default is None.

  • **settings

    etol (float)

    If the energy change between iterations is less than etol, minimization is stopped. Default depends on engine used.

    ftol (float)

    If the magnitude of the global force is less than ftol, minimization is stopped. Default depends on engine used.

    maxiter (int)

    Maximum number of iterations of a single structure relaxation procedure. Default depends on engine used.

    maxeval (int)

    Maximum number of force evaluations to perform. Default depends on engine used.

run(n_steps: int, equilibration: bool = False, verbose: bool = False, output_log: str = None, work_dir: str = None, **settings: dict) None[source]

Runs the MD simulation for the specified number of steps. Trajectories for the simulation are only saved when equilibration is False. Additionally running equilibration for an NVE system (neither barostat nor thermostat set) will temporarily apply a Berendsen thermostat (it is removed from the simulation after the run is completed).

Parameters:
  • n_steps (int) – Number of simulation steps to run

  • equilibration (bool, optional) – If the run is for equilibration (True) or production (False). Default is False.

  • verbose (bool, optional) – Whether to print statements upon starting and completing the run. Default is False.

  • output_log (str, optional) – Log file for the MD engine to write to. Default is None.

  • work_dir (str, optional) – Working directory for the MD engine to write to. Default is None.

property time_step: float

Get or set the simulation time step in fs

Returns:

Simulation time step in fs

Return type:

float

property traj_step: int

Get or set the number of simulation steps between saving the CompactTrajectory

Returns:

Number of simulation steps that elapse between the CompactTrajectory being stored

Return type:

int

property trajectory: CompactTrajectory | None

The CompactTrajectory produced by the most recent production run of the Simulation.

Returns:

Most recent production run CompactTrajectory, or None if no production run has been performed

Return type:

CompactTrajectory

class MDMC.MD.simulation.Universe(dimensions, force_field=None, structures=None, **settings)[source]

Bases: AtomContainer

Class where the configuration and topology are defined

Parameters:
  • dimensions (numpy.ndarray, list, float) – Dimensions of the Universe, in units of Ang. A float can be used for a cubic universe.

  • force_fields (ForceField, optional) – A force field to apply to the Universe. The force fields available are: SPC, OPLSAA, TIP3PFB, SPCE, TIP3P. Default is None.

  • structures (list, optional) – Structure objects contained in the Universe. Default is None.

  • **settings

    kspace_solver (KSpaceSolver)

    The k-space solver to be used for both electrostatic and dispersive interactions. If this is passed then no electrostatic_solver or dispersive_solver may be passed.

    electrostatic_solver (KSpaceSolver)

    The k-space solver to be used for electrostatic interactions.

    dispersive_solver (KSpaceSolver)

    The k-space solver to be used for dispersive interactions.

    constraint_algorithm (ConstraintAlgorithm)

    The constraint algorithm which will be applied to constrained BondedInteractions.

    verbose (str)

    If the output of the class instantiation should be reported, default to True.

dimensions

Dimensions of the Universe in units of Ang.

Type:

numpy.ndarray, list, float

configuration

Stores the content, i.e. configuration of atoms etc within the universe

Type:

Configuration

force_fields

Force field applied to apply to the Universe

Type:

ForceField or None

kspace_solver

The k-space solver to be used for both electrostatic and dispersive interactions.

Type:

KSpaceSolver

electrostatic_solver

The k-space solver to be used for electrostatic interactions.

Type:

KSpaceSolver

dispersive_solver

The k-space solver to be used for dispersive interactions.

Type:

KSpaceSolver

constraint_algorithm

The constraint algorithm which will be applied to constrained BondedInteractions.

Type:

ConstraintAlgorithm

interactions
bonded_interactions
nonbonded_interactions
bonded_interaction_pairs
n_bonded
n_nonbonded
n_interactions
parameters
volume
element_list
element_dict
element_lookup
atoms
n_atoms
molecule_list
n_molecules
structure_list
top_level_structure_list
equivalent_top_level_structures_dict
force_fields
atom_types
atom_type_interactions
density
solvent_density
nbis_by_atom_type_pairs
add_bonded_interaction_pairs(*bonded_interaction_pairs: tuple[Interaction, Atom]) None[source]

Adds one or more interaction pairs to the Universe

Parameters:

*bonded_interaction_pairs – one or more (interaction, atoms) pairs, where atoms is a tuple of all atoms for that specific bonded interaction

add_force_field(force_field: str, *interactions: Interaction, **settings: dict) None[source]

Adds a force field to the specified interactions. If no interactions are passed, the force field is applied to all interactions in the Universe.

Parameters:
  • force_field (str) – The ForceField to parameterize interactions (if provided), or all the interactions in the Universe. The available ForceField are: SPC, OPLSAA, TIP3PFB, SPCE, TIP3P

  • *interactionsInteraction objects to parameterize with the ForceField

  • **settings

    add_dispersions (bool or list of Atoms)

    If True, a Dispersion interaction will be added to all atoms in the Universe. If a list of Atom objects is provided, the Dispersion will be added to these instead. Any added Dispersion interactions (and any previously defined) will then be parametrized by the ForceField. The Dispersion interactions added will only be like-like. By default, no Dispersion interactions are added.

add_nonbonded_interaction(*nonbonded_interactions: tuple[Interaction, Atom]) None[source]

Adds one or more nonbonded interactions to the Universe

Parameters:

*nonbonded_interactions

Nonbonded interactions to be added to the Universe. Can take any number of non-bonded interactions:

Dispersion():

either Lennard-Jones or Buckingham dispersion

Coulombic():

normal or modified Coulomb interaction

with appropriate parameters for the interaction. See http://mdmcproject.org/how-to/use-MDMC/notebooks/defining-molecule-interactions.ipynb for more details on non-bonded interactions.

add_structure(structure: Structure | int, force_field: str = None, center: bool = False) None[source]

Adds a single Structure to the Universe, with optional ForceField applying only to that Structure

Parameters:
  • structure (Structure or int) – The Structure (or its ID as an int) to be added to the Universe

  • force_field (str, optional) – The force field to be applied to the structure. The available ForceField are: SPC, OPLSAA, TIP3PFB, SPCE, TIP3P

  • center (bool, optional) – Whether to center structure within the Universe as it is added

property atom_type_interactions: dict[int, list[Interaction]]

Get the atom types and the interactions for each atom type

Returns:

atom_type:interactions pairs where atom_type is a int specifying the atom type and interactions is a list of Interaction objects acting on that atom_type

Return type:

dict

property atom_types: list[Atom]

Get the atom types of atoms in the Universe

Returns:

The atom types in the Universe

Return type:

list

property atoms: list[Atom]

Get a list of the atoms in the Universe

Returns:

The atoms in the Universe

Return type:

list

property bonded_interaction_pairs: list

Get the bonded interactions and the atoms they apply to

Returns:

The (interaction, atoms) pairs in the Universe, where atoms is a tuple of all Atom objects for that specific interaction

Return type:

list

Example

For an O Atom with two bonds, one to H1 and one to H2:

>>> print(O.bonded_interaction_pairs)
[(Bond, (H1, O)), (Bond, (H2, O))]
property bonded_interactions: list

Get the bonded interactions in the Universe

Returns:

The bonded interactions in the Universe

Return type:

list

property density
property dimensions
property element_dict: dict[str, Atom]

Get the elements in the Universe and example Atom objects for each element

This is required for MD engines which assign the same potential parameters for all identical element types.

Returns:

element:atom pairs, where atom is a single Atom of the specified element.

Return type:

dict

property element_list: list[str]

The elements of the atoms in the Universe

Returns:

The elements in the Universe

Return type:

list

property element_lookup: dict[str, Atom]

Get the elements by ID in the Universe

This is required for MD engines which assign the same potential parameters for all identical element names.

Returns:

element:atom pairs, where atom is a single Atom of the specified element.

Return type:

dict

property equivalent_top_level_structures_dict: dict[Structure, int]

Get a dict of equivalent top level Structure objects that exist in the Universe as keys, and the number of Structure that are equivalent to the key as a value. This does not include any Structure that are a subunit of another structure belonging to the Universe.

Returns:

The top level Structure objects in the Universe as keys, and the number of each as a value

Return type:

dict

fill(structures: Structure, force_field: str = None, num_density: float = None, num_struc_units: int = None) None[source]

A liquid-like filling of the Universe independent of existing atoms

Adds copies of structures to existing configuration until Universe is full. As exclusion region is defined by the size of a bounding sphere, this method is most suitable for atoms or molecules with approximately equal dimensions.

Note

CURRENT APPROACH RESULTS IN NUMBER DENSITY DIFFERENT TO WHAT IS SPECIFIED DEPENDING ON HOW CLOSE CUBE ROOT OF N_MOLECULES IS TO AN int.

Note

CURRENT IMPLEMENTATION SHOULD NOT BE USED WITH NON-CUBIC UNIVERSES AS THE DENSITY MAY OR MAY NOT BE ISOTROPIC DEPENDING ON THE DIMENSIONS AND NUMBER OF UNITS.

Parameters:
  • structures (Structure or int) – The Structure with which to fill the Universe

  • force_field (str) – Applies a ForceField to the Universe. The available ForceField are: SPC, OPLSAA, TIP3PFB, SPCE, TIP3P

  • num_density (float) – Non-negative float specifying the number density of the Structure, in units of Structure / Ang ^ -3

  • num_struc_units (int) – Non-negative int specifying the number of passed Structure objects that the universe should be filled with, regardless of Universe.dimensions.

Raises:
  • ValueError – If both num_density and num_struc_units are passed

  • ValueError – If neither num_density or num_struc_units are passed

property force_fields: None

Get the ForceField acting on the Universe

The available force fields are: SPC, OPLSAA, TIP3PFB, SPCE, TIP3P

Returns:

ForceField that applies to the Universe

Return type:

list

property interactions: list

Get the interactions in the Universe

Returns:

The interactions in the Universe

Return type:

list

property molecule_list: list[Molecule]

Get a list of the Molecule objects in the Universe

Returns:

The Molecule objects in the Universe

Return type:

list

property n_atoms: int

Get the number of atoms in the Universe

Returns:

The number of atoms in the Universe

Return type:

int

property n_bonded: int

Get the number of bonded interactions in the Universe

Returns:

The number of bonded interactions in the Universe

Return type:

int

property n_interactions: int

Get the number of interactions in the Universe

Returns:

The number of interactions in the Universe

Return type:

int

property n_molecules: int

Get the number of molecules in the Universe

Returns:

The number of molecules in the Universe

Return type:

int

property n_nonbonded: int

Get the number of nonbonded interactions in the Universe

Returns:

The number of nonbonded interactions in the Universe

Return type:

int

property nbis_by_atom_type_pairs: dict[tuple[int], list[Interaction]]

Generates a dict of all nonbonded interactions possessed by all combinations of atom_type pairs in the Universe.

Returns:

A dict of {pair: interactions} pairs, where:
  • pair is a tuple for a pair of atom_types in the Universe

  • interactions is a list of nonbonded interactions that exist for this pair of atom_types

Any Dispersions in interactions are ones that exist explicity for this pair, whereas any Coulombics in interactions are ones that exist for either of the atom_types in pair.

Return type:

dict

property nonbonded_interactions: list

Get the nonbonded interactions in the Universe

Returns:

The nonbonded interactions in the Universe

Return type:

list

property parameters: Parameters

Get the parameters of the interactions that exist within the Universe

Returns:

The Parameters objects defined within Universe

Return type:

Parameters

solvate(density: float, tolerance: float = 1.0, solvent: str = 'SPCE', max_iterations: int = 100, **settings: dict) None[source]

Fills the Universe with solvent molecules according to pre-defined coordinates.

Parameters:
  • density (float) – The desired density of the Solvent that solvates the Universe, in units of amu Ang ^ -3

  • tolerance (float, optional) – The +/- percentage tolerance of the density to be achieved. The default is 1 %. Tolerances of less than 1 % are at risk of not converging.

  • solvent (str, optional) – A str specifying an inbuilt Solvent from the following: TIP3PFB, SPCE, SPC, TIP3P. The default is ‘SPCE’.

  • max_iterations (int, optional) – The maximum number of times to try to solvate the universe to within the required density before stopping. Defaults to 100.

  • **settings

    constraint_algorithm (ConstraintAlgorithm)

    A ConstraintAlgorithm which is applied to the Universe. If an inbuilt Solvent is selected (e.g. ‘SPCE’) and constraint_algorithm is not passed, the ConstraintAlgorithm will default to Shake(1e-4, 100).

Raises:

ValueError – If the Universe has already been solvated with a different density.

property solvent_density
property structure_list: None

Get a list of all Structure objects that exist in the Universe. This includes all Structure that are a subunit of another structure belonging to the Universe.

Returns:

The Structure objects in the Universe

Return type:

list

property top_level_structure_list: list[Structure]

Get a list of the top level Structure objects that exist in the Universe. This does not include any Structure that are a subunit of another structure belonging to the Universe.

Returns:

The top level Structure objects in the Universe

Return type:

list

property volume

MDMC.MD.structures module

Module in which all structural units are defined.

Atom is the fundamental structural unit in terms of which all others must be defined. All shared behaviour is included within the Structure base class.

class MDMC.MD.structures.Atom(element: str, position: list[float] | tuple[float] | ndarray = (0.0, 0.0, 0.0), velocity: list[float] | tuple[float] | ndarray = (0.0, 0.0, 0.0), charge: float = None, **settings: dict)[source]

Bases: Structure

A single atom

Parameters:
  • element (periodictable.core.Element) – The periodictable atomic element instance

  • position (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or array with the position in units of Ang. The default is (0., 0., 0.).

  • velocity (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or array with the velocity in units of Ang / fs. Note that if set, the velocity of atoms in the MD engine will be scaled when creating a Simulation in order to ensure the temperature is accurate. Otherwise, if the velocities of all Atom objects in a Simulation are 0, then the velocities of the atoms in the MD engine will be randomly chosen in order to provide an accurate temperature. The default is (0., 0., 0.).

  • charge (float) – The charge of the Atom in units of elementary charge (e). The default is None, meaning that a Coulomb interaction is not applied to the Atom.

  • **settings

    mass (float)

    The atomic mass in amu. If not provided a lookup table will be used.

    atom_type (int)

    All atoms with the same atom_type will have the same non-bonded interactions applied to them. If not provided, MDMC will automatically infer/assign atom types. All atoms with the same element and interactions will have the same atom type.

    name (str)

    A name to uniquely identify the atom. Used by ForceFields (e.g. OPLSAA). Atoms should only have the same names if they are equivalent. Defaults to the element of the atom.

    cutoff (float)

    Sets the cutoff radius in Ang, beyond which this atom does not interact with other atoms. Must be set if a charge is being added to the atom.

element

The periodictable atomic element instance

Type:

periodictable.core.Element

add_interaction(interaction: Interaction, from_interaction: bool = False) None[source]

Adds an interaction to the Atom

Parameters:
  • interaction (MDMC.MD.interactions.Interaction) – Any class dervied from Interaction, or any object with base class Interaction. If an Interaction class is passed then it must be a NonBondedInteraction i.e. only takes a single Atom as an argument. If an Interaction object is passed then this Atom must be in the interaction.atoms.

  • from_interaction (bool, optional) – Specifies if this method has been called from an Interaction. Default is False.

property atom_type: int

Get or set the atom type of the Atom

Returns:

The atom type

Return type:

int

Raises:

AttributeError – The atom_type cannot be changed once it has been set

property atoms: list[Atom]

Get a list of the atoms, just consisting of the Atom

Returns:

A list with a single item, the Atom

Return type:

list

property bonded_interaction_pairs: list

Get bonded interactions acting on the Atom and the other atoms to which the Atom is bonded

Returns:

(interaction, atoms) pairs acting on the Atom, where atoms is a tuple of all Atom objects for that specific bonded interaction.

Return type:

list

Example

For an O Atom with two bonds, one to H1 and one to H2:

>>> print(O.bonded_interaction_pairs)
[(Bond, (H1, O)), (Bond, (H2, O))]
property charge: float | None

Get or set the charge in e if one has been applied to the Atom

If the Atom does not have a Coulombic interaction, setting a value of the charge will create one, and a default cutoff of 10. Ang will be applied

Returns:

The charge in units of e, or None if no charge has been set

Return type:

float

Raises:
  • ValueError – When the Atom has more than one Coulombic interaction

  • ValueError – When the Atom has more than one parameter; i.e. should only have charge as a parameter

  • ValueError – When setting charge to None when a Coulombic interaction already exists.

copy(position: list[float] | tuple[float] | ndarray) Atom[source]

Copies the Atom and all attributes, except ID which is generated

Copying an Atom creates an exact duplicate at the specified position. The copied Atom will have identical bonded and nonbonded interactions as the original. For BondedInteractions this means that the copied atom will be bonded to all atoms to which the original atom is bonded. The ID of the copied atom will differ from the original, as they are sequentially generated.

Parameters:

position (list, tuple, numpy.ndarray) – A 3 element list, tuple or array with the position of the new Atom

Returns:

A copy of the Atom with the specified position

Return type:

Atom

Examples

If the following Atom is copied:

H1 = Atom('H', position=(0., 0., 0.), charge=0.4238)
H2 = H1.copy(position=(1., 1., 1.))

then the new Atom (H2) will have no BondedInteractions, but will have a Coulombic interaction, with a charge of 0.4238 e

If H1 and H2 are then bonded together and a copy is made:

HHbond = Bond((H1, H2))
H3 = H1.copy(position=(2., 2., 2.))

then the newest Atom (H3) will have a Coulombic interaction (also with a charge of 0.4238 e), and it will also have a Bond interaction with H2 (as H1 had a Bond interaction with H2).

copy_interactions(atom: Atom, memo: dict = None) None[source]

This replicates the interactions from self for Atom, but with self substituted by atom in the Interaction.atoms. These interactions are added to any that already exist for the Atom.

Passing the memo dict enables specific interactions to be excluded from being copied, duplicating the behaviour of __deepcopy__

Parameters:
  • atom (Atom) – The Atom for which self.interactions are being replicated

  • memo (dict, optional) – The memoization dict

is_equivalent(structure: Structure) bool[source]

Checks the passed Structure against self for equivalence in terms of the force field application, namely that the following are the same:

  • Element or chemical formula

  • Mass

  • Charge

  • Bonded interactions

  • Non bonded interactions

Returns:

Whether the two are equivalent

Return type:

bool

property mass: float

Get or set the atomic mass in amu

Returns:

the atomic mass in amu

Return type:

float

property nonbonded_interactions: list[NonBondedInteraction]

Get a list of the nonbonded interactions acting on the Atom

Returns:

NonBondedInteractions acting on the Atom

Return type:

list

property universe: Universe | None

Get the Universe to which the Atomm belongs

Returns:

The Universe to which the Atom belongs or None

Return type:

Universe

class MDMC.MD.structures.BoundingBox(atoms: list[Atom])[source]

Bases: object

A box with the minimum and maximum extents of the positions of a collection of atoms

Parameters:

atoms (list) – Atom objects for which the minimum and maximum extents are determined

property max: ndarray

Get or set the maximum extent of the positions of a collection of atoms

Returns:

The maximum extent in Ang

Return type:

numpy.ndarray

property min: ndarray

Get or set the minimum extent of the positions of a collection of atoms

Returns:

The minimum extent in Ang

Return type:

numpy.ndarray

property volume
class MDMC.MD.structures.CompositeStructure(position: list[float] | tuple[float] | ndarray, velocity: list[float] | tuple[float] | ndarray, name: str)[source]

Bases: Structure, AtomContainer

Base class for structural units comprised of more than one Atom

abstract property bonded_interaction_pairs: list[Interaction]

Get bonded interactions acting on the Structure and the other atoms to which the atom is bonded

Returns:

(interaction, atoms) pairs acting on the Structure, where atoms is a tuple of all Atom objects for that specific bonded interaction. At least one of these Atom objects belongs to the Structure

Return type:

list

Example

For an O Atom with two bonds, one to H1 and one to H2:

>>> print(O.bonded_interaction_pairs)
[(Bond, (H1, O)), (Bond, (H2, O))]
calc_CoM() ndarray[source]
Returns:

Position of the center of mass of the CompositeStructure in units of Ang

Return type:

numpy.ndarray

copy(position, rotation=None) CompositeStructure[source]

Copies the CompositeStructure and all attributes, except ID which is generated

Copying a CompositeStructure (e.g. a Molecule) will copy all of the Atom objects within it. All of these atoms will have identical bonded and nonbonded interactions to the CompositeStructure from which they were copied i.e. the CompositeStructure will be exacltly duplicated. The only attributes of the CompositeStructure which will differ are the position (which is passed as a Parameter to copy), and the ID, which is generated automatically.

This will not currently work if the CompositeStructure has any bonded interactions with atoms external to it (e.g. it may cause issues for copying molecules with groups)

Interactions for Atom objects may be reordered with respect to initial atoms

Parameters:
  • position (list, tuple, numpy.ndarray) – 3 element list, tuple or array of float specifying the position of the new Structure

  • rotation (list, tuple, numpy.ndarray, optional) – 3 element list, tuple or array of floats specifying the degrees of anticlockwise rotation around the x, y, and z axes respectively. The rotation is centered on the center of mass of the CompositeStructure. The default rotation is None, which applies no rotation to the copied CompositeStructure.

Returns:

A CompositeStructure of the same type with all non-unique attributes copied and a new position

Return type:

CompositeStructure

property formula: str

Get the chemical formula of the CompositeStructure

Returns:

The chemical formula using the Hill system

Return type:

str

abstract property nonbonded_interactions: list[Interaction]

Get a list of the nonbonded interactions acting on the Structure

Returns:

NonBondedInteractions acting on the Structure

Return type:

list

rotate(x: float = 0.0, y: float = 0.0, z: float = 0.0) None[source]

Rotates the CompositeStructure around its center of mass

In all cases (e.g. x, y and z) the rotation is anticlockwise about the specific axis

Parameters:
  • x (float, optional) – The angle of rotation around the x-axis in deg. The default is 0.

  • y (float, optional) – The angle of rotation around the y-axis in deg. The default is 0.

  • z (float, optional) – The angle of rotation around the z-axis in deg. The default is 0.

property structure_list: list[Structure]

Get or set the Structure objects that are subunits of this CompositeStructure

Returns:

list of Structure that are subunits of this CompositeStructure

Return type:

list

property universe: Universe | None

Get or set the Universe to which the CompositeStructure belongs

Returns:

The Universe to which the CompositeStructure belongs or None

Return type:

Universe

class MDMC.MD.structures.Molecule(position: list[float] | tuple[float] | ndarray = None, velocity: list[float] | tuple[float] | ndarray = (0, 0, 0), name=None, **settings: dict)[source]

Bases: CompositeStructure

Two or more bonded atoms

Must be declared with at least 2 Atom objects.

Parameters:
  • position (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or array with the position in units of Ang. The default is None, which sets the position of the Molecule to be equal to the center of mass of the atoms in the Molecule.

  • velocity (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or array with the velocity in units of Ang / fs. The default is (0., 0., 0.).

  • name (str, optional) – The name of the structure. The default is None.

  • **settings

    atoms (list)

    A list of Atom which will be included in the Molecule

    interactions (list)

    A list of Interaction acting on atoms within the Molecule. The interactions provides a convenience for declaring interactions on atoms when a Molecule is initialized. It is not required and is exactly equivalent to initializing the interactions prior to the Molecule.

property bonded_interaction_pairs: list

Get bonded interactions acting on the Molecule

Returns:

(interaction, atoms) pairs acting on the Molecule, where atoms is a tuple of all atoms for that specific bonded interaction. At least one of these atoms belongs to the Molecule

Return type:

list

Example

For an O Atom with two bonds, one to H1 and one to H2:

>>> print(O.bonded_interaction_pairs)
[(Bond, (H1, O)), (Bond, (H2, O))]
property charge
is_equivalent(structure: Structure) bool[source]

Checks the passed Structure against self for equivalence in terms of the force field application, namely that the following are the same:

  • Element or chemical formula

  • Mass

  • Charge

  • Bonded interactions

  • Non bonded interactions

Returns:

Whether the two are equivalent

Return type:

bool

property mass
property nonbonded_interactions: list[NonBondedInteraction]

Get a list of the nonbonded interactions acting on the Molecule

Returns:

NonBondedInteraction objects acting on the Molecule

Return type:

list

property position: ndarray

Get or set the position of the center of mass of the Molecule in Ang

Also sets the positions of all subunits

Return type:

numpy.ndarray

class MDMC.MD.structures.Structure(position: list[float] | tuple[float] | ndarray, velocity: list[float] | tuple[float] | ndarray, name: str)[source]

Bases: ABC

Abstract base class for all structural units

Parameters:
  • position (list, tuple, numpy.ndarray) – A 3 element list, tuple or array with the position in units of Ang.

  • velocity (list, tuple, numpy.ndarray) – A 3 element list, tuple or array with the velocity in units of Ang / fs.

  • name (str) –

    The name of the structure.

    Attributes

  • ---------- – ID : int A unique identifier for each Structure.

  • universe (Universe) – The Universe to which the Structure belongs.

  • name – The name of the structure.

  • parent (Structure) – Structure to which this unit belongs, or self

property atoms: list[Atom]

Get a list of all of the Atom objects in the structure by recursively calling atoms for all substructures

Returns:

All atoms in the structure

Return type:

list

abstract property bonded_interaction_pairs: list[tuple[Interaction, tuple[Atom]]]

Get bonded interactions acting on the Structure and the other atoms to which the atom is bonded

Returns:

(interaction, atoms) pairs acting on the Structure, where atoms is a tuple of all Atom objects for that specific bonded interaction. At least one of these Atom objects belongs to the Structure

Return type:

list

Example

For an O Atom with two bonds, one to H1 and one to H2:

>>> print(O.bonded_interaction_pairs)
[(Bond, (H1, O)), (Bond, (H2, O))]
property bonded_interactions: list[BondedInteraction]

Get a list of the bonded interactions acting on the Structure

Returns:

BondedInteractions acting on the Structure

Return type:

list

property bounding_box: BoundingBox

returns: Contains the lower and upper extents of the Molecule :rtype: BoundingBox

copy(position: list[float] | tuple[float] | ndarray) Structure[source]

Copies the structural unit and sets the position

Parameters:

position (list, tuple, numpy.ndarray) – 3 element list, tuple or array of float specifying the position of the new Structure

Returns:

A Structure of the same type with all non-unique attributes copied and a new position

Return type:

Structure

property interactions: list[Interaction]

Get a list of the interactions acting on the Structure

Returns:

Interactions acting on the Structure

Return type:

list

abstract is_equivalent(structure: Structure) bool[source]

Checks the passed Structure against self for equivalence in terms of the force field application, namely that the following are the same:

  • Element or chemical formula

  • Mass

  • Charge

  • Bonded interactions

  • Non bonded interactions

Returns:

Whether the two are equivalent

Return type:

bool

abstract property nonbonded_interactions: list[NonBondedInteraction]

Get a list of the nonbonded interactions acting on the Structure

Returns:

NonBondedInteractions acting on the Structure

Return type:

list

property position: ndarray

Get or set the position of the center of mass of the Structure in Ang

Return type:

numpy.ndarray

property structure_type: str

Get the class of the Structure.

Returns:

The name of the class

Return type:

str

property top_level_structure: Structure

Get the top level structure (i.e. Structure which has no parent) of the Structure

Returns:

Highest level Structure of which it is a member

Return type:

Structure

translate(displacement: tuple | ndarray) None[source]

Translate the structural unit by the specified displacement

Parameters:

Displacement (tuple, numpy.ndarray) – A three element tuple or array of float

abstract property universe: 'Universe' | None

Get the Universe to which the Structure belongs

Returns:

The Universe to which the Structure belongs or None

Return type:

Universe

valid_position(position: list[float] | tuple[float] | ndarray = None) bool[source]

Checks if the specified position is within the bounds of the Structure.universe, if it has one

Parameters:

position (list, tuple, numpy.ndarray) – 3 element list, tuple or array with units of Ang or None. If None then the position of the Structure is used.

Returns:

True if position is within Universe or there is no associated Universe. False if Structure has an associated Universe but the position is not within its bounds.

Return type:

bool

Raises:

ValueError – If position if undefined

property velocity: ndarray

Get or set the velocity of the Structure in Ang/fs

Return type:

numpy.ndarray

MDMC.MD.structures.filter_atoms(atoms: list[Atom], predicate: Callable) list[Atom][source]

Filters a list of Atoms with a given predicate

Parameters:
  • atoms (list) – A list of Atom

  • predicate (function) – A function that returns a bool

Returns:

Atom objects in atoms which meet the condition of predicate

Return type:

list

MDMC.MD.structures.filter_atoms_element(atoms: list[Atom], element: str) list[Atom][source]

Filters a list of atoms based on the atomic element

Parameters:
  • atoms (list) – A list of Atom

  • element (str) – The atomic element label

Returns:

Atom objects of a specific element

Return type:

list

MDMC.MD.structures.get_reduced_chemical_formula(symbols: list[str], factor: int = None, system: str = 'Hill') str[source]

Get the reduced chemical formula

Parameters:
  • symbols (list of str) – The chemical formula to be reduced. It is expressed as a list of elements, with a single element for each atom. Elements are grouped by type but not ordered e.g. all 'O' values, then all 'H' values etc.

  • factor (int, optional) – The factor by which the total number of symbols will be reduced. If None, the greatest common divisor of the different symbols will be used. The default is None.

  • system (str, optional) – Determines the order of the chemical formula. If 'Hill' the Hill system is used to determine the order. If None, the order is based on the order of `symbols. The default is 'Hill'.

Returns:

The chemical formula corresponding to symbols, except with only n_atoms. If system is 'Hill', the formula will be ordered as per the Hill system, otherwise the formula will be ordered based on the order of symbols.

Return type:

str

Example

Reducing the formula for four water molecules to a single water molecules:

>>> get_reduced_chemical_formula(['H'] * 8 + ['O'] * 4)
'H2O'

Reducing the formula for four water molecules to two water molecules:

>>> get_reduced_chemical_formula(['H'] * 8 + ['O'] * 4, factor=2)
'H4O2'

Module contents

Classes and functions for setting up and running a molecular dynamics simulation

Contents

ase engine_facades force_fields interaction_functions solvents simulation structures interactions