MDMC.MD package
Subpackages
- MDMC.MD.ase package
- MDMC.MD.engine_facades package
- Submodules
- MDMC.MD.engine_facades.dlpoly_engine module
DLPOLYAttributeDLPOLYEngineDLPOLYEngine.HANDLED_PARAMSDLPOLYEngine.barostatDLPOLYEngine.clear()DLPOLYEngine.convert_trajectory()DLPOLYEngine.ensembleDLPOLYEngine.eval()DLPOLYEngine.minimize()DLPOLYEngine.pressureDLPOLYEngine.reset_config()DLPOLYEngine.run()DLPOLYEngine.save_config()DLPOLYEngine.saved_configDLPOLYEngine.setup_simulation()DLPOLYEngine.setup_universe()DLPOLYEngine.temperatureDLPOLYEngine.thermostatDLPOLYEngine.update_parameters()
DLPOLYEnsembleDLPOLYSimulationDLPOLYUniverseconvert_unit()
- MDMC.MD.engine_facades.facade module
MDEngineMDEngine.clear()MDEngine.convert_trajectory()MDEngine.eval()MDEngine.minimize()MDEngine.parent_simulationMDEngine.reset_config()MDEngine.run()MDEngine.save_config()MDEngine.saved_configMDEngine.setup_simulation()MDEngine.setup_universe()MDEngine.time_stepMDEngine.traj_stepMDEngine.update_parameters()
MDEngineError
- MDMC.MD.engine_facades.facade_factory module
- MDMC.MD.engine_facades.lammps_engine module
LAMMPSEngineLAMMPSEngine.barostatLAMMPSEngine.clear()LAMMPSEngine.convert_trajectory()LAMMPSEngine.ensembleLAMMPSEngine.eval()LAMMPSEngine.minimize()LAMMPSEngine.pressureLAMMPSEngine.reset_config()LAMMPSEngine.run()LAMMPSEngine.save_config()LAMMPSEngine.saved_configLAMMPSEngine.setup_simulation()LAMMPSEngine.setup_universe()LAMMPSEngine.temperatureLAMMPSEngine.thermostatLAMMPSEngine.update_parameters()
LAMMPSEnsembleLAMMPSEnsemble.rescale_stepLAMMPSEnsemble.apply_ensemble_fixes()LAMMPSEnsemble.barostatLAMMPSEnsemble.p_dampLAMMPSEnsemble.pressureLAMMPSEnsemble.remove_ensemble_fixes()LAMMPSEnsemble.t_dampLAMMPSEnsemble.t_fractionLAMMPSEnsemble.t_windowLAMMPSEnsemble.temperatureLAMMPSEnsemble.thermostatLAMMPSEnsemble.time_step
LAMMPSSimulationLAMMPSSimulation.universeLAMMPSSimulation.traj_stepLAMMPSSimulation.ensembleLAMMPSSimulation.ang_momentum_stepsLAMMPSSimulation.barostatLAMMPSSimulation.lin_momentum_stepsLAMMPSSimulation.neighbor_stepsLAMMPSSimulation.pressureLAMMPSSimulation.skinLAMMPSSimulation.temperatureLAMMPSSimulation.thermostatLAMMPSSimulation.time_step
LAMMPSUniverseLAMMPSUniverse.universeLAMMPSUniverse.atom_dictLAMMPSUniverse.atom_typesLAMMPSUniverse.atom_type_propertiesLAMMPSUniverse.bondsLAMMPSUniverse.anglesLAMMPSUniverse.coulsLAMMPSUniverse.dispsLAMMPSUniverse.bond_IDLAMMPSUniverse.angle_IDLAMMPSUniverse.proper_IDLAMMPSUniverse.improper_IDLAMMPSUniverse.apply_constraints()LAMMPSUniverse.nonbonded_mixLAMMPSUniverse.set_config()LAMMPSUniverse.update_parameters()
PyLammpsAttributeconvert_unit()parse_all_nonbonded_styles()parse_bonded_coefficients()parse_bonded_styles()parse_constraint()parse_dispersion_coefficients()parse_kspace_solver()parse_nonbonded_modifications()parse_nonbonded_styles()
- Module contents
- MDMC.MD.force_fields package
- Submodules
- MDMC.MD.force_fields.OPLSAA module
- MDMC.MD.force_fields.SPC module
- MDMC.MD.force_fields.SPCE module
- MDMC.MD.force_fields.TIP3P module
- MDMC.MD.force_fields.TIP3PFB module
- MDMC.MD.force_fields.ff module
- MDMC.MD.force_fields.force_field_factory module
- Module contents
- MDMC.MD.packmol package
- Submodules
- MDMC.MD.packmol.packmol_setup module
PackmolSetupPackmolSetup.add_box()PackmolSetup.add_container()PackmolSetup.add_cube()PackmolSetup.add_fixed_structure()PackmolSetup.add_sphere()PackmolSetup.get_max_sizes()PackmolSetup.get_settings()PackmolSetup.get_structures()PackmolSetup.remove_structure()PackmolSetup.resolve_density()PackmolSetup.tolerancePackmolSetup.validate_setup()
calculate_volume()
- MDMC.MD.packmol.packmol_wrapper module
- Module contents
PackmolFillerPackmolSetupPackmolSetup.add_box()PackmolSetup.add_container()PackmolSetup.add_cube()PackmolSetup.add_fixed_structure()PackmolSetup.add_sphere()PackmolSetup.get_max_sizes()PackmolSetup.get_settings()PackmolSetup.get_structures()PackmolSetup.remove_structure()PackmolSetup.resolve_density()PackmolSetup.tolerancePackmolSetup.validate_setup()
- MDMC.MD.solvents package
Submodules
MDMC.MD.container module
Module for the AtomContainer class
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:
InteractionFunctionThe Buckingham potential (in units of
kJ mol^-1) for the interaction of 2 atoms at distance r (inAng) 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:
Examples
The following creates a
Dispersioninteraction with aBuckinghamfunctional 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:
InteractionFunctionCoulomb 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
Coulombicinteraction with aCoulombfunctional form:O = Atom('O') coul = Coulomb(-0.8476) O_coulombic = Coulombic(atoms=O, cutoff=10., function=coul)
As
Coulombis the default functional form ofCoulombicinteractions when anAtomis created, this is equivalent to setting thechargeon anAtom:O = Atom('O', charge=-0.8476)
- class MDMC.MD.interaction_functions.HarmonicPotential(equilibrium_state: float, potential_strength: float, **settings: dict)[source]
Bases:
InteractionFunctionHarmonic potential for bond stretching, and angular and improper dihedral vibration, with the form:
\[E = K(r-r_0)^2\]As
HarmonicPotentialcan be used with several differentInteractiontypes, theInteractiontype must be specified, so that the correct units can be assigned to theequilibrium_stateandpotential_strengthparameters.- Parameters:
equilibrium_state (float) – The equilibrium state of the object in either
Angordegrees, depending on theinteraction_typepassed.potential_strength (float) – The potential strength in units of
kJ mol^-1 Ang^-2(linear) orkJ mol^-1 rad^-2(angular/improper), depending on theinteraction_typepassed.**settings –
- interaction_typestr
A str specifying either
'bond','angle'or'improper'. This assigns the correct units toequilibrium_stateandpotential_strength. This keyword must be passed.
- Raises:
ValueError – The
interaction_typemust be'bond','angle', or'improper'TypeError – An
interaction_typeof'bond','angle', or'improper'must be passed
Examples
The following creates a
HarmonicPotentialfor aBondinteraction:hp = HarmonicPotential(1., 4637., interaction_type='bond')
The following creates a
HarmonicPotentialfor aBondAngleinteraction:hp = HarmonicPotential(109.47, 383., interaction_type='angle')
The following creates a
HarmonicPotentialfor aDihedralAngleinteraction withimproper==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:
objectBase class for interaction functions, which can be user supplied
- Parameters:
val_dict (dict) –
name:valuepairs. Currently this must be ordered alphabetically.valuemust either be an object with avalueand aunit(e.g. aUnitFloatobject), 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:
- property parameters: Parameters
Get or set the interaction function’s parameters
- Returns:
A Parameters object containing each
Parameter- Return type:
- property parameters_values: ndarray
Get the values for all
Parametersobjects- Returns:
A NumPy
arrayof values for allParameter- Return type:
- set_parameters_interactions(interaction: Interaction) None[source]
Sets the
parentInteractionfor allParametersobjects- Parameters:
interaction (Interaction) – An
Interactionto set as theparentof all theParameters
- class MDMC.MD.interaction_functions.LennardJones(epsilon: float, sigma: float, **settings: dict)[source]
Bases:
InteractionFunctionDispersive 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^-1sigma (float) – The LJ sigma value in units of
Ang**settings –
- cutofffloat
The distance in
Angat 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
Dispersioninteraction with aLennardJonesfunctional 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:
InteractionFunctionPeriodic 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^-1n1 (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
Periodicfor aDihedralAngleinteraction withimproper==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
Periodicfor aDihedralAngleinteraction withimproper==False(i.e. a proper dihedral), withK1=3.53548,K2=-4.02501andK3=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 ofInteractionFunction.- 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 (orUnit) 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 inparameter_units.
Examples
The following adds units of
'Ang'to parameteralpha, units of's'to the parameterbeta, and units of'atm'to the parametergamma:@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
UnitFloatorUnitNDArraywill not be created). So to setepsilonas 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,BondedInteractionA bond between any two atoms. Requires exactly two atoms in each
atom_tuple.- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomwhich 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 samefunctiondescribing the interaction, and have the sameconstrainedsetting.- Parameters:
other (BondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.interactions.BondAngle(*atom_tuples: tuple, **settings: dict)[source]
Bases:
ConstrainableMixin,BondedInteractionA bond angle between any two bonds
Requires three
Atomobjects (rotation around central atom) in eachatom_tuple. The atoms are orderedi,j,k, wherejis the central atom. So:BondAngle(i, j, k) == BondAngle(k, j, i)
- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomwhich are bonded together. For three or moreAtomobjects, the order of theAtomobjects 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 samefunctiondescribing the interaction, and have the sameconstrainedsetting.- Parameters:
other (BondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.interactions.BondedInteraction(*atom_tuples: list[tuple], **settings: dict)[source]
Bases:
InteractionBase class for bonded interactions
- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomobjects which are bonded together. For three or moreAtomobjects, the order of theAtomobjects within each tuple is important.**settings –
n_atoms(int)The number of atoms to which this
BondedInteractionapplies, for example 2 for aBond.
Examples
For a single bonded interactions which applies to
H1,O1, andH2:BondedInteraction(H1, O1, H2)
For two bonded interactions of the same
BondedInteractiontype, one applied toH1,O1andH2Atomobjects, and the other applied toH3,O2andH4Atomobjects: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-OBondedInteraction: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
Atomobjects**settings –
from_structure(bool)If
add_atomshas been called from aStructure
- Raises:
ValueError – If this
BondedInteractionhas already been applied to one or more of theatoms
- property atom_types: Set[int | None]
Get the set of all
atom_type``s that this ``BondedInteractioncorresponds to, including None if appropriate.
- property atoms: list[tuple[Atom]]
Get or set the atoms on which the
Coulombicinteraction is applied
- element_list() list | None[source]
Get a list of the elements for which the
BondedInteractionapplies- Returns:
The elements for which the
BondedInteractionapplies- Return type:
- is_equivalent(other) bool[source]
Checks for equivalence between two
BondedInteraction``s, specifically if they apply to the same ``atom_typesand the samefunctiondescribing the interaction.- Parameters:
other (BondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.interactions.ConstrainableMixin(*atom_tuples: tuple, **settings: dict)[source]
Bases:
objectA mixin class enabling classes inheriting from
BondedInteractionto be constrainedThese constraints are then applied by a constraint algorithm (e.g. SHAKE), which is specified in the
Universeto which theBondedInteractionbelongs.- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomobjects which are bonded together. For three or moreAtomobjects, the order of theAtomobjects within each tuple is important.**settings
- class MDMC.MD.interactions.Coulombic(universe: Universe = None, **settings: dict)[source]
Bases:
NonBondedInteractionA non-bonded coulombic interaction - either normal or modified Coulomb
- Parameters:
universe (Universe, optional) – The
Universein which theCoulombicexists. Default is None. Must be passed as a parameter ifatom_typesif passed.**settings –
charge(float)The charge parameter of the
Coulombicinteraction, in units ofe. If this argument is passed, theinteraction_functionof thisCoulombicis set toCoulombwith this float as itsParameter. Passingchargewill overwrite any otherinteraction_functionsthat are set, i.e. it makesfunctionparameter redundantatoms(list)Atomobjects to which theCoulombicapplies. If specifying theatoms,universedoes 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:
Examples
Upon initializing an
Atomobject and adding it to aUniverse: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
atomsis passed then aUniversedoes not need to be passed:O_coulombic = Coulombic(atoms=[O], charge=-0.84)
- property atom_types: list
Get the atom types for which the
Coulombicapplies- Returns:
All atom types to which the
Coulombicapplies. If the interaction was initialized withatoms, allatom_typesof theatomsto which theCoulombicwas applied are returned; HOWEVER THE COULOMBIC INTERACTION IS NOT APPLIED TO ALL ATOMS OF THESEatom_types, ONLY THE ATOMS INself.atoms- Return type:
- class MDMC.MD.interactions.DihedralAngle(*atom_tuples: tuple, **settings: dict)[source]
Bases:
BondedInteractionA dihedral angle between any two sets of three atoms,
ijkandjkl.Dihedral angles can be both proper and improper, where the angle between the two planes of
ijkandjklis fixed for improper dihedrals.The atoms of a proper
DihedralAngleare orderedi,j,k,l, wherejandkare the two central atoms. So:DihedralAngle(i, j, k, l) == DihedralAngle(l, k, j, i)
The atoms of an improper
DihedralAngleare orderedi,j,k,l, whereiis the central atom to whichj,k, andlare 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
DihedralAngleis improper or not.
- class MDMC.MD.interactions.Dispersion(universe: Universe, *atom_types: int, **settings: dict)[source]
Bases:
NonBondedInteractionA non-bonded dispersive interaction - either LJ or Buckingham
- Parameters:
universe (Universe) – The
Universein which theNonBondedInteractionexists*atom_types – int for each atom type for which the
NonBondedInteractionapplies**settings –
cutoff(float)The distance in
Angat which the interaction potential is truncatedvdw_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:
TypeError –
atom_typesmust be iterableValueError –
Dispersionshould only be specified as existing between pairs ofatom_typesTypeError – Each
atom_typemust be int
- property atom_types
Get the atom types for which the
NonBondedInteractionapplies- Returns:
A list of int for the
atom_types- Return type:
- property atoms: list[tuple[list[Atom]]]
Get the atoms on which the
Dispersionis applied- Returns:
A list of two tuple, where each tuple contains a list of Atom. Every
Atomin the first tuple has a dispersion interaction with everyAtomin 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:
- element_list() list[source]
Get a list of the elements for which the
Interactionapplies- Returns:
The elements for which the
Interactionapplies- Return type:
- is_equivalent(other) bool[source]
Checks for equivalence between two
Dispersion``s, specifically if they apply to the same ``atom_types, have the samecuttoff, the samefunctiondescribing the interaction andvdw_tail_correctionsetting.- Parameters:
other (Dispersion) – The object to compare against.
- Return type:
- class MDMC.MD.interactions.Interaction(**settings: dict)[source]
Bases:
ABCBase class for interactions, both bonded, non-bonded and constraints
Each different type of interaction should have an
Interactionobject. This object contains a list of theAtom(or`Atompairs, triplets or quadruplets, depending on the type of interaction) for which thisInteractionapplies. For example, an oxygenCoulombicinteraction would contain a list of tuple where each tuple contains a different OAtom, and a hydrogen-oxygenBondinteraction would contain a list of tuple where each tuple contains a different H and O pair.Interactionobjects can be sliced to return a sublist of the tuple.When an
Atomis passed to anInteraction, theInteractionis also added to theAtom.- Parameters:
**settings –
function(InteractionFunction)A class of interaction function (e.g.
HarmonicPotential)
- abstract element_list() list[source]
Get a list of the elements for which the
Interactionapplies- Returns:
The elements for which the
Interactionapplies- Return type:
- element_tuple() tuple[source]
A tuple of elements for which the
Interactionapplies- Returns:
The elements for which the
Interactionapplies- Return type:
- property function: InteractionFunction
Get or set the
InteractionFunctionof theInteraction- Returns:
The interaction function of the
Interaction- Return type:
- property function_name: str | None
Get the name of the
InteractionFunctionbelonging to theInteraction- Returns:
The name of the
InteractionFunction, or None if noInteractionFunctionhas been set- Return type:
- property parameters: Parameters
Get the
Parameterobjects belonging to theInteractionFunctionbelonging to theInteraction- Returns:
A
Parametersobject containing eachParameter- Return type:
- class MDMC.MD.interactions.NonBondedInteraction(universe, *atom_types, **settings)[source]
Bases:
InteractionBase class for non-bonded interactions
- Parameters:
universe (Universe) – The
Universein which theNonBondedInteractionexists*atom_types – int for each
atom_typefor which theNonBondedInteractionapplies**settings –
cutoff(float)The distance in
Angat which the interaction potential is truncated
- abstract property atom_types: list[int]
Get the atom types for which the
NonBondedInteractionapplies- Returns:
A list of int for the
atom_types- Return type:
- property cutoff: float
Get or set the distance in
Angat which the interaction potential is truncated- Returns:
The distance in
Angof thecutoff- Return type:
- is_equivalent(other) bool[source]
Checks for equivalence between two
NonBondedInteraction``s, specifically if they apply to the same ``atom_types, have the samecuttoffand the samefunctiondescribing the interaction.- Parameters:
other (NonBondedInteraction) – The object to compare against.
- Return type:
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:
objectA 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).constraintsmust have the same units asvalue.**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
Interactionobjects for thisParameter- Returns:
All parent
Interactionobjects- Return type:
- Raises:
ValueError – If an added interaction name is not consistent with existing interaction names
ValueError – If an added
Interactionhas a function name not consistent with the function names of an existingInteraction
- set_tie(parameter: Parameter, expr: str) None[source]
This
tiestheParameter.valueto thevalueof anotherParameterExamples
To set the
Parameter.valuetop1.value * 2:>>> Parameter.set_tie(p1, "* 2")
- property tie: float | None
Get the
valueof a theParameterthat thisParameteris tied to- Returns:
The
valueof thetiedParameter- Return type:
- property tied: bool
Get whether this
Parameteris tied- Returns:
True if this
Parameteris tied to anotherParameter, else False- Return type:
- static validate_value(value: float, constraints: tuple) None[source]
Validates the
Parameter.valueby testing if it is within theconstraints- Parameters:
- Raises:
ValueError – If the
valueis not within theconstraints
- property value: float
Get or set the value of the
ParameterThe value will not be changed if it is
fixedortied, or if it is set outside the bounds ofconstraints- Returns:
The value of the
Parameter, including if theParameteristied- Return type:
- Warns:
warnings.warn – If the
Parameterisfixed.warnings.warn – If the
Parameteristied.
- class MDMC.MD.parameters.Parameters(init_parameters: Parameter | list[Parameter] | None = None)[source]
Bases:
dictA dict-like object where every element is a
Parameterindexed by name, which contains a number of helper methods for filtering.Although
Parametersis 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 (
Parameteror list of ``Parameter``s, optional, default None) – The initialParameterobjects that theParametersobject 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
Parameteror list of ``Parameter``s to the dict, with the parameter name as its key.- Parameters:
parameters (
Parameteror 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
Parameteras an argument.- Returns:
The
Parameterobjects which meet the condition of the predicate- Return type:
- filter_atom_attribute(attribute: str, value: str | float) Parameters[source]
Filters based on the attribute of
Atomobjects which have eachParameterapplied to them- Parameters:
- Returns:
The
Parameterobjects which are applied to anAtomobject which has the specifiedvalueof the specifiedattribute- Return type:
- filter_function(function_name: str) Parameters[source]
Filters based on the name of the
InteractionFunctionof eachParameter- Parameters:
function_name (str) – The name of the
InteractionFunctionofParameterobjects to return, for example'LennardJones'or'HarmonicPotential'.- Returns:
The
Parameterobjects which have afunctionwith the specifiedfunction_name- Return type:
- filter_interaction(interaction_name: str) Parameters[source]
Filters based on the name of the
Interactionof eachParameter- Parameters:
interaction_name (str) – The name of the
InteractionofParameterobjects to return, for example'Bond'.- Returns:
The
Parameterobjects which have anInteractionwith the specifiedinteraction_name- Return type:
- filter_name(name: str) Parameters[source]
Filters by
name- Parameters:
name (str) – The
nameof theParameterobjects to return.- Returns:
The
Parameterobjects withname- Return type:
- filter_structure(structure_name: str) Parameters[source]
Filters based on the name of the
Structureto which eachParameterapplies- Parameters:
structure_name (str) – The name of a
Structure.- Returns:
The
Parameterobjects which are applied to aStructurewhich has the specifiedstructure_name- Return type:
- filter_value(comparison: str, value: float) Parameters[source]
Filters by
value- Parameters:
- Returns:
The
Parameterobjects which return a True when their values are compared withvalueusing thecomparisonoperator- Return type:
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:
objectClass describing the algorithm and parameters which are applied to constrain
BondedInteractionobjects- Parameters:
- class MDMC.MD.simulation.Ewald(**settings: dict)[source]
Bases:
KSpaceSolverHolds 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:
objectClass describing the k-space solver that is applied to electrostatic and/or dispersion interactions
Different
MDEnginerequire 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
- class MDMC.MD.simulation.PPPM(**settings: dict)[source]
Bases:
KSpaceSolverHolds 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:
ConstraintAlgorithmHolds the parameters which are required for the RATTLE algorithm to be applied to the constrained interactions
- class MDMC.MD.simulation.Shake(accuracy: float, max_iterations: int)[source]
Bases:
ConstraintAlgorithmHolds the parameters which are required for the SHAKE algorithm to be applied to the constrained interactions
- class MDMC.MD.simulation.Simulation(universe: Universe, traj_step: int, time_step: float = 1.0, engine: str = 'lammps', **settings)[source]
Bases:
objectSets 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
Universeon which the simulation is performed.traj_step (int) – How many steps the simulation should take between dumping each
CompactTrajectoryframe. Along withtime_stepdetermines the time separation of calculated variables such as energy.time_step (float, optional) – Simulation timestep in
fs. Default is 1.engine (str, optional) – The
MDEngineused 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 thetemperature. If all atoms in theuniversehave 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 theengine, not theuniverse, are affected.integrator(str)Simulation time integrator.
lj_options(dict)option:valuepairs for Lennard-Jones interactions.es_options(dict)option:valuepairs 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.
- traj_step
How many steps the simulation should take between dumping each
CompactTrajectoryframe. Along withtime_stepdetermines the time separation of calculated variables such as energy.- Type:
- engine
A subclass of
MDEnginewhich provides the interface to the MD library. Default is'lammps'.- Type:
MDEngine, optional
- time_step
- trajectory
- auto_equilibrate(variables: list[str] = ('temp', 'pe'), eq_step: int = 10, window_size: int = 100, tolerance: float = 0.01, max_step=10000) 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
equilibrationis 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
CompactTrajectorybeing stored- Return type:
int
- property trajectory: CompactTrajectory | None
The
CompactTrajectoryproduced by the most recent production run of theSimulation.- Returns:
Most recent production run
CompactTrajectory, or None if no production run has been performed- Return type:
- class MDMC.MD.simulation.Universe(dimensions, force_field=None, structures=None, **settings)[source]
Bases:
AtomContainerClass where the configuration and topology are defined
- Parameters:
dimensions (numpy.ndarray, list, float) – Dimensions of the
Universe, in units ofAng. 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: OPLSAA, TIP3P, TIP3PFB, SPCE, SPC. Default is None.
structures (list, optional) –
Structureobjects contained in theUniverse. 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_solverordispersive_solvermay 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
Universein units ofAng.- Type:
- configuration
Stores the content, i.e. configuration of atoms etc within the universe
- Type:
- 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:
- electrostatic_solver
The k-space solver to be used for electrostatic interactions.
- Type:
- dispersive_solver
The k-space solver to be used for dispersive interactions.
- Type:
- constraint_algorithm
The constraint algorithm which will be applied to constrained
BondedInteractions.- Type:
- 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, whereatomsis a tuple of all atoms for that specific bondedinteraction
- add_force_field(force_field: str, *interactions: Interaction, **settings: dict) None[source]
Adds a force field to the specified
interactions. If nointeractionsare passed, the force field is applied to all interactions in theUniverse.- Parameters:
force_field (str) – The
ForceFieldto parameterizeinteractions(if provided), or all theinteractionsin theUniverse. The availableForceFieldare: OPLSAA, TIP3P, TIP3PFB, SPCE, SPC*interactions –
Interactionobjects to parameterize with theForceField**settings –
add_dispersions(bool or list ofAtoms)If True, a
Dispersioninteraction will be added to all atoms in theUniverse. If a list ofAtomobjects is provided, theDispersionwill be added to these instead. Any addedDispersioninteractions (and any previously defined) will then be parametrized by theForceField. TheDispersioninteractions added will only be like-like. By default, noDispersioninteractions 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
Structureto theUniverse, with optionalForceFieldapplying only to thatStructure- Parameters:
structure (Structure or int) – The
Structure(or itsIDas an int) to be added to theUniverseforce_field (str, optional) – The force field to be applied to the structure. The available
ForceFieldare: OPLSAA, TIP3P, TIP3PFB, SPCE, SPCcenter (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:interactionspairs whereatom_typeis a int specifying the atom type andinteractionsis a list ofInteractionobjects acting on thatatom_type- Return type:
- property atom_types: list[Atom]
Get the atom types of atoms in the
Universe- Returns:
The atom types in the
Universe- Return type:
- property atoms: list[Atom]
Get a list of the atoms in the
Universe- Returns:
The atoms in the
Universe- Return type:
- property bonded_interaction_pairs: list
Get the bonded interactions and the atoms they apply to
- Returns:
The (
interaction,atoms) pairs in theUniverse, whereatomsis a tuple of allAtomobjects for that specificinteraction- Return type:
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:
- property density
- property dimensions
- property element_dict: dict[str, Atom]
Get the elements in the
Universeand exampleAtomobjects for each elementThis is required for MD engines which assign the same potential parameters for all identical element types.
- Returns:
element:atom pairs, whereatomis a singleAtomof the specifiedelement.- Return type:
- property element_list: list[str]
The elements of the atoms in the
Universe- Returns:
The elements in the
Universe- Return type:
- property element_lookup: dict[str, Atom]
Get the elements by ID in the
UniverseThis is required for MD engines which assign the same potential parameters for all identical element names.
- Returns:
element:atom pairs, whereatomis a singleAtomof the specifiedelement.- Return type:
- property equivalent_top_level_structures_dict: dict[Structure, int]
Get a dict of equivalent top level
Structureobjects that exist in theUniverseas keys, and the number ofStructurethat are equivalent to the key as a value. This does not include anyStructurethat are a subunit of another structure belonging to theUniverse.- Returns:
The top level
Structureobjects in theUniverseas keys, and the number of each as a value- Return type:
- fill(structures: Structure, force_field: str = None, num_density: float = None, num_struc_units: int = None) None[source]
A liquid-like filling of the
Universeindependent of existing atomsAdds copies of
structuresto existing configuration untilUniverseis 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
Structurewith which to fill theUniverseforce_field (str) – Applies a
ForceFieldto theUniverse. The availableForceFieldare: OPLSAA, TIP3P, TIP3PFB, SPCE, SPCnum_density (float) – Non-negative float specifying the number density of the
Structure, in units ofStructure / Ang ^ -3num_struc_units (int) – Non-negative int specifying the number of passed
Structureobjects that the universe should be filled with, regardless ofUniverse.dimensions.
- Raises:
ValueError – If both
num_densityandnum_struc_unitsare passedValueError – If neither
num_densityornum_struc_unitsare passed
- property force_fields: None
Get the
ForceFieldacting on theUniverseThe available force fields are: OPLSAA, TIP3P, TIP3PFB, SPCE, SPC
- Returns:
ForceFieldthat applies to theUniverse- Return type:
- property interactions: list
Get the interactions in the
Universe- Returns:
The interactions in the
Universe- Return type:
- property molecule_list: list[Molecule]
Get a list of the
Moleculeobjects in theUniverse- Returns:
The
Moleculeobjects in theUniverse- Return type:
- property n_atoms: int
Get the number of atoms in the
Universe- Returns:
The number of atoms in the
Universe- Return type:
- property n_bonded: int
Get the number of bonded interactions in the
Universe- Returns:
The number of bonded interactions in the
Universe- Return type:
- property n_interactions: int
Get the number of interactions in the
Universe- Returns:
The number of interactions in the
Universe- Return type:
- property n_molecules: int
Get the number of molecules in the
Universe- Returns:
The number of molecules in the
Universe- Return type:
- property n_nonbonded: int
Get the number of nonbonded interactions in the
Universe- Returns:
The number of nonbonded interactions in the
Universe- Return type:
- property nbis_by_atom_type_pairs: dict[tuple[int], list[Interaction]]
Generates a dict of all nonbonded interactions possessed by all combinations of
atom_typepairs in theUniverse.- Returns:
- A dict of
{pair: interactions}pairs, where: pairis a tuple for a pair ofatom_typesin theUniverseinteractionsis a list of nonbonded interactions that exist for thispairofatom_types
Any
Dispersionsininteractionsare ones that exist explicity for thispair, whereas anyCoulombicsininteractionsare ones that exist for either of theatom_typesinpair.- A dict of
- Return type:
- property nonbonded_interactions: list
Get the nonbonded interactions in the
Universe- Returns:
The nonbonded interactions in the
Universe- Return type:
- property parameters: Parameters
Get the parameters of the interactions that exist within the
Universe- Returns:
The
Parametersobjects defined withinUniverse- Return type:
- solvate(density: float, tolerance: float = 1.0, solvent: str = 'SPCE', max_iterations: int = 100, **settings: dict) None[source]
Fills the
Universewith solvent molecules according to pre-defined coordinates.- Parameters:
density (float) – The desired density of the
Solventthat solvates theUniverse, in units ofamu Ang ^ -3tolerance (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
Solventfrom 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
ConstraintAlgorithmwhich is applied to theUniverse. If an inbuiltSolventis selected (e.g. ‘SPCE’) andconstraint_algorithmis not passed, theConstraintAlgorithmwill default toShake(1e-4, 100).
- Raises:
ValueError – If the
Universehas already been solvated with a different density.
- property solvent_density
- property structure_list: None
Get a list of all
Structureobjects that exist in theUniverse. This includes allStructurethat are a subunit of another structure belonging to theUniverse.- Returns:
The
Structureobjects in theUniverse- Return type:
- property top_level_structure_list: list[Structure]
Get a list of the top level
Structureobjects that exist in theUniverse. This does not include anyStructurethat are a subunit of another structure belonging to theUniverse.- Returns:
The top level
Structureobjects in theUniverse- Return type:
- 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:
StructureA single atom
- Parameters:
element (periodictable.core.Element) – The periodictable atomic element instance
position (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the position in units ofAng. The default is(0., 0., 0.).velocity (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the velocity in units ofAng / 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
Atomin units of elementary charge (e). The default is None, meaning that aCoulombinteraction is not applied to theAtom.**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
- 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 classInteraction. If anInteractionclass is passed then it must be aNonBondedInteractioni.e. only takes a singleAtomas an argument. If anInteractionobject is passed then thisAtommust be in theinteraction.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:
- Raises:
AttributeError – The
atom_typecannot 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:
- property bonded_interaction_pairs: list
Get bonded interactions acting on the
Atomand the other atoms to which theAtomis bonded- Returns:
(
interaction,atoms) pairs acting on theAtom, whereatomsis a tuple of all Atom objects for that specific bondedinteraction.- Return type:
Example
For an O
Atomwith 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
eif one has been applied to theAtomIf the
Atomdoes not have aCoulombicinteraction, setting a value of thechargewill create one, and a defaultcutoffof10. Angwill be applied- Returns:
The charge in units of
e, or None if no charge has been set- Return type:
- Raises:
ValueError – When the
Atomhas more than oneCoulombicinteractionValueError – When the
Atomhas more than one parameter; i.e. should only have charge as a parameterValueError – When setting charge to None when a
Coulombicinteraction already exists.
- copy(position: list[float] | tuple[float] | ndarray) Atom[source]
Copies the
Atomand all attributes, exceptIDwhich is generatedCopying an
Atomcreates an exact duplicate at the specifiedposition. The copiedAtomwill have identical bonded and nonbonded interactions as the original. ForBondedInteractionsthis means that the copied atom will be bonded to all atoms to which the original atom is bonded. TheIDof 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
arraywith thepositionof the newAtom- Returns:
A copy of the
Atomwith the specifiedposition- Return type:
Examples
If the following
Atomis 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 noBondedInteractions, but will have aCoulombicinteraction, with achargeof0.4238 eIf
H1andH2are 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 aCoulombicinteraction (also with achargeof0.4238 e), and it will also have aBondinteraction withH2(asH1had aBondinteraction withH2).
- copy_interactions(atom: Atom, memo: dict = None) None[source]
This replicates the interactions from
selfforAtom, but withselfsubstituted byatomin theInteraction.atoms. These interactions are added to any that already exist for theAtom.Passing the
memodict enables specific interactions to be excluded from being copied, duplicating the behaviour of__deepcopy__
- is_equivalent(structure: Structure) bool[source]
Checks the passed
Structureagainst 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:
- property mass: float
Get or set the atomic mass in
amu- Returns:
the atomic mass in
amu- Return type:
- property nonbonded_interactions: list[NonBondedInteraction]
Get a list of the nonbonded interactions acting on the
Atom- Returns:
NonBondedInteractionsacting on theAtom- Return type:
- class MDMC.MD.structures.BoundingBox(atoms: list[Atom])[source]
Bases:
objectA box with the minimum and maximum extents of the positions of a collection of atoms
- Parameters:
atoms (list) –
Atomobjects 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:
- 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:
- property volume
- class MDMC.MD.structures.CompositeStructure(position: list[float] | tuple[float] | ndarray, velocity: list[float] | tuple[float] | ndarray, name: str)[source]
Bases:
Structure,AtomContainerBase class for structural units comprised of more than one
Atom- abstract property bonded_interaction_pairs: list[Interaction]
Get bonded interactions acting on the
Structureand the other atoms to which the atom is bonded- Returns:
(
interaction,atoms) pairs acting on theStructure, whereatomsis a tuple of allAtomobjects for that specific bondedinteraction. At least one of theseAtomobjects belongs to theStructure- Return type:
Example
For an O
Atomwith 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
CompositeStructurein units ofAng- Return type:
- copy(position, rotation=None) CompositeStructure[source]
Copies the
CompositeStructureand all attributes, exceptIDwhich is generatedCopying a
CompositeStructure(e.g. aMolecule) will copy all of theAtomobjects within it. All of these atoms will have identical bonded and nonbonded interactions to theCompositeStructurefrom which they were copied i.e. theCompositeStructurewill be exacltly duplicated. The only attributes of theCompositeStructurewhich will differ are theposition(which is passed as a Parameter tocopy), and theID, which is generated automatically.This will not currently work if the
CompositeStructurehas any bonded interactions with atoms external to it (e.g. it may cause issues for copying molecules with groups)Interactions for
Atomobjects may be reordered with respect to initial atoms- Parameters:
position (list, tuple, numpy.ndarray) – 3 element list, tuple or
arrayof float specifying thepositionof the newStructurerotation (list, tuple, numpy.ndarray, optional) – 3 element list, tuple or
arrayof 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 theCompositeStructure. The defaultrotationis None, which applies no rotation to the copiedCompositeStructure.
- Returns:
A
CompositeStructureof the same type with all non-unique attributes copied and a newposition- Return type:
- property formula: str
Get the chemical formula of the
CompositeStructure- Returns:
The chemical formula using the Hill system
- Return type:
- abstract property nonbonded_interactions: list[Interaction]
Get a list of the nonbonded interactions acting on the
Structure- Returns:
NonBondedInteractionsacting on theStructure- Return type:
- rotate(x: float = 0.0, y: float = 0.0, z: float = 0.0) None[source]
Rotates the
CompositeStructurearound its center of massIn all cases (e.g. x, y and z) the rotation is anticlockwise about the specific axis
- 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:
CompositeStructureTwo or more bonded atoms
Must be declared with at least 2
Atomobjects.- Parameters:
position (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the position in units ofAng. The default is None, which sets the position of theMoleculeto be equal to the center of mass of the atoms in theMolecule.velocity (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the velocity in units ofAng / 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
Atomwhich will be included in theMoleculeinteractions(list)A list of
Interactionacting on atoms within theMolecule. Theinteractionsprovides a convenience for declaring interactions on atoms when aMoleculeis initialized. It is not required and is exactly equivalent to initializing the interactions prior to theMolecule.
- property bonded_interaction_pairs: list
Get bonded interactions acting on the
Molecule- Returns:
(
interaction,atoms) pairs acting on theMolecule, whereatomsis a tuple of all atoms for that specific bondedinteraction. At least one of theseatomsbelongs to theMolecule- Return type:
Example
For an
OAtomwith two bonds, one toH1and one toH2:>>> print(O.bonded_interaction_pairs) [(Bond, (H1, O)), (Bond, (H2, O))]
- property charge
- is_equivalent(structure: Structure) bool[source]
Checks the passed
Structureagainst 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:
- property mass
- property nonbonded_interactions: list[NonBondedInteraction]
Get a list of the nonbonded interactions acting on the
Molecule- Returns:
NonBondedInteractionobjects acting on theMolecule- Return type:
- class MDMC.MD.structures.Structure(position: list[float] | tuple[float] | ndarray, velocity: list[float] | tuple[float] | ndarray, name: str)[source]
Bases:
ABCAbstract base class for all structural units
- Parameters:
position (list, tuple, numpy.ndarray) – A 3 element list, tuple or
arraywith the position in units ofAng.velocity (list, tuple, numpy.ndarray) – A 3 element list, tuple or
arraywith the velocity in units ofAng / fs.name (str) –
The name of the structure.
Attributes
---------- – ID : int A unique identifier for each
Structure.universe (Universe) – The
Universeto which theStructurebelongs.name – The name of the structure.
parent (Structure) –
Structureto which this unit belongs, orself
- property atoms: list[Atom]
Get a list of all of the Atom objects in the structure by recursively calling
atomsfor all substructures- Returns:
All atoms in the structure
- Return type:
- abstract property bonded_interaction_pairs: list[tuple[Interaction, tuple[Atom]]]
Get bonded interactions acting on the
Structureand the other atoms to which the atom is bonded- Returns:
(
interaction,atoms) pairs acting on theStructure, whereatomsis a tuple of allAtomobjects for that specific bondedinteraction. At least one of theseAtomobjects belongs to theStructure- Return type:
Example
For an O
Atomwith 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:
BondedInteractionsacting on theStructure- Return type:
- 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
arrayof float specifying thepositionof the newStructure- Returns:
A
Structureof the same type with all non-unique attributes copied and a newposition- Return type:
- property interactions: list[Interaction]
Get a list of the interactions acting on the
Structure- Returns:
Interactions acting on the
Structure- Return type:
- abstract is_equivalent(structure: Structure) bool[source]
Checks the passed
Structureagainst 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:
- abstract property nonbonded_interactions: list[NonBondedInteraction]
Get a list of the nonbonded interactions acting on the
Structure- Returns:
NonBondedInteractionsacting on theStructure- Return type:
- property position: ndarray
Get or set the position of the center of mass of the
StructureinAng- Return type:
- property structure_type: str
Get the class of the
Structure.- Returns:
The name of the class
- Return type:
- property top_level_structure: Structure
Get the top level structure (i.e.
Structurewhich has noparent) of theStructure- Returns:
Highest level
Structureof which it is a member- Return type:
- translate(displacement: tuple | ndarray) None[source]
Translate the structural unit by the specified displacement
- Parameters:
Displacement (tuple, numpy.ndarray) – A three element tuple or
arrayof float
- abstract property universe: 'Universe' | None
Get the
Universeto which theStructurebelongs- Returns:
The
Universeto which theStructurebelongs or None- Return type:
- valid_position(position: list[float] | tuple[float] | ndarray = None) bool[source]
Checks if the specified
positionis within the bounds of theStructure.universe, if it has one- Parameters:
position (list, tuple, numpy.ndarray) – 3 element list, tuple or
arraywith units ofAngor None. If None then thepositionof theStructureis used.- Returns:
True if
positionis withinUniverseor there is no associatedUniverse. False ifStructurehas an associatedUniversebut thepositionis not within its bounds.- Return type:
- Raises:
ValueError – If
positionif undefined
- MDMC.MD.structures.filter_atoms(atoms: list[Atom], predicate: Callable) list[Atom][source]
Filters a list of Atoms with a given predicate
- MDMC.MD.structures.filter_atoms_element(atoms: list[Atom], element: str) list[Atom][source]
Filters a list of atoms based on the atomic element
- 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 onlyn_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 ofsymbols.- Return type:
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
- class MDMC.MD.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:
StructureA single atom
- Parameters:
element (periodictable.core.Element) – The periodictable atomic element instance
position (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the position in units ofAng. The default is(0., 0., 0.).velocity (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the velocity in units ofAng / 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
Atomin units of elementary charge (e). The default is None, meaning that aCoulombinteraction is not applied to theAtom.**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
- 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 classInteraction. If anInteractionclass is passed then it must be aNonBondedInteractioni.e. only takes a singleAtomas an argument. If anInteractionobject is passed then thisAtommust be in theinteraction.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:
- Raises:
AttributeError – The
atom_typecannot 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:
- property bonded_interaction_pairs: list
Get bonded interactions acting on the
Atomand the other atoms to which theAtomis bonded- Returns:
(
interaction,atoms) pairs acting on theAtom, whereatomsis a tuple of all Atom objects for that specific bondedinteraction.- Return type:
Example
For an O
Atomwith 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
eif one has been applied to theAtomIf the
Atomdoes not have aCoulombicinteraction, setting a value of thechargewill create one, and a defaultcutoffof10. Angwill be applied- Returns:
The charge in units of
e, or None if no charge has been set- Return type:
- Raises:
ValueError – When the
Atomhas more than oneCoulombicinteractionValueError – When the
Atomhas more than one parameter; i.e. should only have charge as a parameterValueError – When setting charge to None when a
Coulombicinteraction already exists.
- copy(position: list[float] | tuple[float] | ndarray) Atom[source]
Copies the
Atomand all attributes, exceptIDwhich is generatedCopying an
Atomcreates an exact duplicate at the specifiedposition. The copiedAtomwill have identical bonded and nonbonded interactions as the original. ForBondedInteractionsthis means that the copied atom will be bonded to all atoms to which the original atom is bonded. TheIDof 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
arraywith thepositionof the newAtom- Returns:
A copy of the
Atomwith the specifiedposition- Return type:
Examples
If the following
Atomis 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 noBondedInteractions, but will have aCoulombicinteraction, with achargeof0.4238 eIf
H1andH2are 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 aCoulombicinteraction (also with achargeof0.4238 e), and it will also have aBondinteraction withH2(asH1had aBondinteraction withH2).
- copy_interactions(atom: Atom, memo: dict = None) None[source]
This replicates the interactions from
selfforAtom, but withselfsubstituted byatomin theInteraction.atoms. These interactions are added to any that already exist for theAtom.Passing the
memodict enables specific interactions to be excluded from being copied, duplicating the behaviour of__deepcopy__
- is_equivalent(structure: Structure) bool[source]
Checks the passed
Structureagainst 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:
- property mass: float
Get or set the atomic mass in
amu- Returns:
the atomic mass in
amu- Return type:
- property nonbonded_interactions: list[NonBondedInteraction]
Get a list of the nonbonded interactions acting on the
Atom- Returns:
NonBondedInteractionsacting on theAtom- Return type:
- class MDMC.MD.Bond(*atom_tuples: tuple, **settings: dict)[source]
Bases:
ConstrainableMixin,BondedInteractionA bond between any two atoms. Requires exactly two atoms in each
atom_tuple.- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomwhich 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 samefunctiondescribing the interaction, and have the sameconstrainedsetting.- Parameters:
other (BondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.BondAngle(*atom_tuples: tuple, **settings: dict)[source]
Bases:
ConstrainableMixin,BondedInteractionA bond angle between any two bonds
Requires three
Atomobjects (rotation around central atom) in eachatom_tuple. The atoms are orderedi,j,k, wherejis the central atom. So:BondAngle(i, j, k) == BondAngle(k, j, i)
- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomwhich are bonded together. For three or moreAtomobjects, the order of theAtomobjects 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 samefunctiondescribing the interaction, and have the sameconstrainedsetting.- Parameters:
other (BondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.BondedInteraction(*atom_tuples: list[tuple], **settings: dict)[source]
Bases:
InteractionBase class for bonded interactions
- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomobjects which are bonded together. For three or moreAtomobjects, the order of theAtomobjects within each tuple is important.**settings –
n_atoms(int)The number of atoms to which this
BondedInteractionapplies, for example 2 for aBond.
Examples
For a single bonded interactions which applies to
H1,O1, andH2:BondedInteraction(H1, O1, H2)
For two bonded interactions of the same
BondedInteractiontype, one applied toH1,O1andH2Atomobjects, and the other applied toH3,O2andH4Atomobjects: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-OBondedInteraction: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
Atomobjects**settings –
from_structure(bool)If
add_atomshas been called from aStructure
- Raises:
ValueError – If this
BondedInteractionhas already been applied to one or more of theatoms
- property atom_types: Set[int | None]
Get the set of all
atom_type``s that this ``BondedInteractioncorresponds to, including None if appropriate.
- property atoms: list[tuple[Atom]]
Get or set the atoms on which the
Coulombicinteraction is applied
- element_list() list | None[source]
Get a list of the elements for which the
BondedInteractionapplies- Returns:
The elements for which the
BondedInteractionapplies- Return type:
- is_equivalent(other) bool[source]
Checks for equivalence between two
BondedInteraction``s, specifically if they apply to the same ``atom_typesand the samefunctiondescribing the interaction.- Parameters:
other (BondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.BoundingBox(atoms: list[Atom])[source]
Bases:
objectA box with the minimum and maximum extents of the positions of a collection of atoms
- Parameters:
atoms (list) –
Atomobjects 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:
- 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:
- property volume
- class MDMC.MD.Buckingham(A: float, B: float, C: float)[source]
Bases:
InteractionFunctionThe Buckingham potential (in units of
kJ mol^-1) for the interaction of 2 atoms at distance r (inAng) 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:
Examples
The following creates a
Dispersioninteraction with aBuckinghamfunctional 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.CompositeStructure(position: list[float] | tuple[float] | ndarray, velocity: list[float] | tuple[float] | ndarray, name: str)[source]
Bases:
Structure,AtomContainerBase class for structural units comprised of more than one
Atom- abstract property bonded_interaction_pairs: list[Interaction]
Get bonded interactions acting on the
Structureand the other atoms to which the atom is bonded- Returns:
(
interaction,atoms) pairs acting on theStructure, whereatomsis a tuple of allAtomobjects for that specific bondedinteraction. At least one of theseAtomobjects belongs to theStructure- Return type:
Example
For an O
Atomwith 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
CompositeStructurein units ofAng- Return type:
- copy(position, rotation=None) CompositeStructure[source]
Copies the
CompositeStructureand all attributes, exceptIDwhich is generatedCopying a
CompositeStructure(e.g. aMolecule) will copy all of theAtomobjects within it. All of these atoms will have identical bonded and nonbonded interactions to theCompositeStructurefrom which they were copied i.e. theCompositeStructurewill be exacltly duplicated. The only attributes of theCompositeStructurewhich will differ are theposition(which is passed as a Parameter tocopy), and theID, which is generated automatically.This will not currently work if the
CompositeStructurehas any bonded interactions with atoms external to it (e.g. it may cause issues for copying molecules with groups)Interactions for
Atomobjects may be reordered with respect to initial atoms- Parameters:
position (list, tuple, numpy.ndarray) – 3 element list, tuple or
arrayof float specifying thepositionof the newStructurerotation (list, tuple, numpy.ndarray, optional) – 3 element list, tuple or
arrayof 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 theCompositeStructure. The defaultrotationis None, which applies no rotation to the copiedCompositeStructure.
- Returns:
A
CompositeStructureof the same type with all non-unique attributes copied and a newposition- Return type:
- property formula: str
Get the chemical formula of the
CompositeStructure- Returns:
The chemical formula using the Hill system
- Return type:
- abstract property nonbonded_interactions: list[Interaction]
Get a list of the nonbonded interactions acting on the
Structure- Returns:
NonBondedInteractionsacting on theStructure- Return type:
- rotate(x: float = 0.0, y: float = 0.0, z: float = 0.0) None[source]
Rotates the
CompositeStructurearound its center of massIn all cases (e.g. x, y and z) the rotation is anticlockwise about the specific axis
- class MDMC.MD.ConstrainableMixin(*atom_tuples: tuple, **settings: dict)[source]
Bases:
objectA mixin class enabling classes inheriting from
BondedInteractionto be constrainedThese constraints are then applied by a constraint algorithm (e.g. SHAKE), which is specified in the
Universeto which theBondedInteractionbelongs.- Parameters:
atom_tuples (list) – A list of tuple. Each tuple contains
Atomobjects which are bonded together. For three or moreAtomobjects, the order of theAtomobjects within each tuple is important.**settings
- class MDMC.MD.ConstraintAlgorithm(accuracy: float, max_iterations: int)[source]
Bases:
objectClass describing the algorithm and parameters which are applied to constrain
BondedInteractionobjects- Parameters:
- class MDMC.MD.Coulomb(charge: float)[source]
Bases:
InteractionFunctionCoulomb 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
Coulombicinteraction with aCoulombfunctional form:O = Atom('O') coul = Coulomb(-0.8476) O_coulombic = Coulombic(atoms=O, cutoff=10., function=coul)
As
Coulombis the default functional form ofCoulombicinteractions when anAtomis created, this is equivalent to setting thechargeon anAtom:O = Atom('O', charge=-0.8476)
- class MDMC.MD.Coulombic(universe: Universe = None, **settings: dict)[source]
Bases:
NonBondedInteractionA non-bonded coulombic interaction - either normal or modified Coulomb
- Parameters:
universe (Universe, optional) – The
Universein which theCoulombicexists. Default is None. Must be passed as a parameter ifatom_typesif passed.**settings –
charge(float)The charge parameter of the
Coulombicinteraction, in units ofe. If this argument is passed, theinteraction_functionof thisCoulombicis set toCoulombwith this float as itsParameter. Passingchargewill overwrite any otherinteraction_functionsthat are set, i.e. it makesfunctionparameter redundantatoms(list)Atomobjects to which theCoulombicapplies. If specifying theatoms,universedoes 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:
Examples
Upon initializing an
Atomobject and adding it to aUniverse: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
atomsis passed then aUniversedoes not need to be passed:O_coulombic = Coulombic(atoms=[O], charge=-0.84)
- property atom_types: list
Get the atom types for which the
Coulombicapplies- Returns:
All atom types to which the
Coulombicapplies. If the interaction was initialized withatoms, allatom_typesof theatomsto which theCoulombicwas applied are returned; HOWEVER THE COULOMBIC INTERACTION IS NOT APPLIED TO ALL ATOMS OF THESEatom_types, ONLY THE ATOMS INself.atoms- Return type:
- class MDMC.MD.DihedralAngle(*atom_tuples: tuple, **settings: dict)[source]
Bases:
BondedInteractionA dihedral angle between any two sets of three atoms,
ijkandjkl.Dihedral angles can be both proper and improper, where the angle between the two planes of
ijkandjklis fixed for improper dihedrals.The atoms of a proper
DihedralAngleare orderedi,j,k,l, wherejandkare the two central atoms. So:DihedralAngle(i, j, k, l) == DihedralAngle(l, k, j, i)
The atoms of an improper
DihedralAngleare orderedi,j,k,l, whereiis the central atom to whichj,k, andlare 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
DihedralAngleis improper or not.
- class MDMC.MD.Dispersion(universe: Universe, *atom_types: int, **settings: dict)[source]
Bases:
NonBondedInteractionA non-bonded dispersive interaction - either LJ or Buckingham
- Parameters:
universe (Universe) – The
Universein which theNonBondedInteractionexists*atom_types – int for each atom type for which the
NonBondedInteractionapplies**settings –
cutoff(float)The distance in
Angat which the interaction potential is truncatedvdw_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:
TypeError –
atom_typesmust be iterableValueError –
Dispersionshould only be specified as existing between pairs ofatom_typesTypeError – Each
atom_typemust be int
- property atom_types
Get the atom types for which the
NonBondedInteractionapplies- Returns:
A list of int for the
atom_types- Return type:
- property atoms: list[tuple[list[Atom]]]
Get the atoms on which the
Dispersionis applied- Returns:
A list of two tuple, where each tuple contains a list of Atom. Every
Atomin the first tuple has a dispersion interaction with everyAtomin 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:
- element_list() list[source]
Get a list of the elements for which the
Interactionapplies- Returns:
The elements for which the
Interactionapplies- Return type:
- is_equivalent(other) bool[source]
Checks for equivalence between two
Dispersion``s, specifically if they apply to the same ``atom_types, have the samecuttoff, the samefunctiondescribing the interaction andvdw_tail_correctionsetting.- Parameters:
other (Dispersion) – The object to compare against.
- Return type:
- class MDMC.MD.Ewald(**settings: dict)[source]
Bases:
KSpaceSolverHolds 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.HarmonicPotential(equilibrium_state: float, potential_strength: float, **settings: dict)[source]
Bases:
InteractionFunctionHarmonic potential for bond stretching, and angular and improper dihedral vibration, with the form:
\[E = K(r-r_0)^2\]As
HarmonicPotentialcan be used with several differentInteractiontypes, theInteractiontype must be specified, so that the correct units can be assigned to theequilibrium_stateandpotential_strengthparameters.- Parameters:
equilibrium_state (float) – The equilibrium state of the object in either
Angordegrees, depending on theinteraction_typepassed.potential_strength (float) – The potential strength in units of
kJ mol^-1 Ang^-2(linear) orkJ mol^-1 rad^-2(angular/improper), depending on theinteraction_typepassed.**settings –
- interaction_typestr
A str specifying either
'bond','angle'or'improper'. This assigns the correct units toequilibrium_stateandpotential_strength. This keyword must be passed.
- Raises:
ValueError – The
interaction_typemust be'bond','angle', or'improper'TypeError – An
interaction_typeof'bond','angle', or'improper'must be passed
Examples
The following creates a
HarmonicPotentialfor aBondinteraction:hp = HarmonicPotential(1., 4637., interaction_type='bond')
The following creates a
HarmonicPotentialfor aBondAngleinteraction:hp = HarmonicPotential(109.47, 383., interaction_type='angle')
The following creates a
HarmonicPotentialfor aDihedralAngleinteraction withimproper==True(i.e. an improper dihedral):hp = HarmonicPotential(180., 20.92, interaction_type='improper')
- class MDMC.MD.Interaction(**settings: dict)[source]
Bases:
ABCBase class for interactions, both bonded, non-bonded and constraints
Each different type of interaction should have an
Interactionobject. This object contains a list of theAtom(or`Atompairs, triplets or quadruplets, depending on the type of interaction) for which thisInteractionapplies. For example, an oxygenCoulombicinteraction would contain a list of tuple where each tuple contains a different OAtom, and a hydrogen-oxygenBondinteraction would contain a list of tuple where each tuple contains a different H and O pair.Interactionobjects can be sliced to return a sublist of the tuple.When an
Atomis passed to anInteraction, theInteractionis also added to theAtom.- Parameters:
**settings –
function(InteractionFunction)A class of interaction function (e.g.
HarmonicPotential)
- abstract element_list() list[source]
Get a list of the elements for which the
Interactionapplies- Returns:
The elements for which the
Interactionapplies- Return type:
- element_tuple() tuple[source]
A tuple of elements for which the
Interactionapplies- Returns:
The elements for which the
Interactionapplies- Return type:
- property function: InteractionFunction
Get or set the
InteractionFunctionof theInteraction- Returns:
The interaction function of the
Interaction- Return type:
- property function_name: str | None
Get the name of the
InteractionFunctionbelonging to theInteraction- Returns:
The name of the
InteractionFunction, or None if noInteractionFunctionhas been set- Return type:
- property parameters: Parameters
Get the
Parameterobjects belonging to theInteractionFunctionbelonging to theInteraction- Returns:
A
Parametersobject containing eachParameter- Return type:
- class MDMC.MD.InteractionFunction(val_dict: dict)[source]
Bases:
objectBase class for interaction functions, which can be user supplied
- Parameters:
val_dict (dict) –
name:valuepairs. Currently this must be ordered alphabetically.valuemust either be an object with avalueand aunit(e.g. aUnitFloatobject), 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:
- property parameters: Parameters
Get or set the interaction function’s parameters
- Returns:
A Parameters object containing each
Parameter- Return type:
- property parameters_values: ndarray
Get the values for all
Parametersobjects- Returns:
A NumPy
arrayof values for allParameter- Return type:
- set_parameters_interactions(interaction: Interaction) None[source]
Sets the
parentInteractionfor allParametersobjects- Parameters:
interaction (Interaction) – An
Interactionto set as theparentof all theParameters
- class MDMC.MD.KSpaceSolver(**settings: dict)[source]
Bases:
objectClass describing the k-space solver that is applied to electrostatic and/or dispersion interactions
Different
MDEnginerequire 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
- class MDMC.MD.LennardJones(epsilon: float, sigma: float, **settings: dict)[source]
Bases:
InteractionFunctionDispersive 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^-1sigma (float) – The LJ sigma value in units of
Ang**settings –
- cutofffloat
The distance in
Angat 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
Dispersioninteraction with aLennardJonesfunctional form:lj = LennardJones(0.6502, 3.166) O_disp = Disperion(universe, (O.atom_type, O.atom_type), function=lj)
- class MDMC.MD.Molecule(position: list[float] | tuple[float] | ndarray = None, velocity: list[float] | tuple[float] | ndarray = (0, 0, 0), name=None, **settings: dict)[source]
Bases:
CompositeStructureTwo or more bonded atoms
Must be declared with at least 2
Atomobjects.- Parameters:
position (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the position in units ofAng. The default is None, which sets the position of theMoleculeto be equal to the center of mass of the atoms in theMolecule.velocity (list, tuple, numpy.ndarray, optional) – A 3 element list, tuple or
arraywith the velocity in units ofAng / 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
Atomwhich will be included in theMoleculeinteractions(list)A list of
Interactionacting on atoms within theMolecule. Theinteractionsprovides a convenience for declaring interactions on atoms when aMoleculeis initialized. It is not required and is exactly equivalent to initializing the interactions prior to theMolecule.
- property bonded_interaction_pairs: list
Get bonded interactions acting on the
Molecule- Returns:
(
interaction,atoms) pairs acting on theMolecule, whereatomsis a tuple of all atoms for that specific bondedinteraction. At least one of theseatomsbelongs to theMolecule- Return type:
Example
For an
OAtomwith two bonds, one toH1and one toH2:>>> print(O.bonded_interaction_pairs) [(Bond, (H1, O)), (Bond, (H2, O))]
- property charge
- is_equivalent(structure: Structure) bool[source]
Checks the passed
Structureagainst 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:
- property mass
- property nonbonded_interactions: list[NonBondedInteraction]
Get a list of the nonbonded interactions acting on the
Molecule- Returns:
NonBondedInteractionobjects acting on theMolecule- Return type:
- class MDMC.MD.NonBondedInteraction(universe, *atom_types, **settings)[source]
Bases:
InteractionBase class for non-bonded interactions
- Parameters:
universe (Universe) – The
Universein which theNonBondedInteractionexists*atom_types – int for each
atom_typefor which theNonBondedInteractionapplies**settings –
cutoff(float)The distance in
Angat which the interaction potential is truncated
- abstract property atom_types: list[int]
Get the atom types for which the
NonBondedInteractionapplies- Returns:
A list of int for the
atom_types- Return type:
- property cutoff: float
Get or set the distance in
Angat which the interaction potential is truncated- Returns:
The distance in
Angof thecutoff- Return type:
- is_equivalent(other) bool[source]
Checks for equivalence between two
NonBondedInteraction``s, specifically if they apply to the same ``atom_types, have the samecuttoffand the samefunctiondescribing the interaction.- Parameters:
other (NonBondedInteraction) – The object to compare against.
- Return type:
- class MDMC.MD.PPPM(**settings: dict)[source]
Bases:
KSpaceSolverHolds 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.Parameter(value, name, fixed=False, constraints=None, **settings)[source]
Bases:
objectA 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).constraintsmust have the same units asvalue.**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
Interactionobjects for thisParameter- Returns:
All parent
Interactionobjects- Return type:
- Raises:
ValueError – If an added interaction name is not consistent with existing interaction names
ValueError – If an added
Interactionhas a function name not consistent with the function names of an existingInteraction
- set_tie(parameter: Parameter, expr: str) None[source]
This
tiestheParameter.valueto thevalueof anotherParameterExamples
To set the
Parameter.valuetop1.value * 2:>>> Parameter.set_tie(p1, "* 2")
- property tie: float | None
Get the
valueof a theParameterthat thisParameteris tied to- Returns:
The
valueof thetiedParameter- Return type:
- property tied: bool
Get whether this
Parameteris tied- Returns:
True if this
Parameteris tied to anotherParameter, else False- Return type:
- static validate_value(value: float, constraints: tuple) None[source]
Validates the
Parameter.valueby testing if it is within theconstraints- Parameters:
- Raises:
ValueError – If the
valueis not within theconstraints
- property value: float
Get or set the value of the
ParameterThe value will not be changed if it is
fixedortied, or if it is set outside the bounds ofconstraints- Returns:
The value of the
Parameter, including if theParameteristied- Return type:
- Warns:
warnings.warn – If the
Parameterisfixed.warnings.warn – If the
Parameteristied.
- class MDMC.MD.Parameters(init_parameters: Parameter | list[Parameter] | None = None)[source]
Bases:
dictA dict-like object where every element is a
Parameterindexed by name, which contains a number of helper methods for filtering.Although
Parametersis 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 (
Parameteror list of ``Parameter``s, optional, default None) – The initialParameterobjects that theParametersobject 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
Parameteror list of ``Parameter``s to the dict, with the parameter name as its key.- Parameters:
parameters (
Parameteror 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
Parameteras an argument.- Returns:
The
Parameterobjects which meet the condition of the predicate- Return type:
- filter_atom_attribute(attribute: str, value: str | float) Parameters[source]
Filters based on the attribute of
Atomobjects which have eachParameterapplied to them- Parameters:
- Returns:
The
Parameterobjects which are applied to anAtomobject which has the specifiedvalueof the specifiedattribute- Return type:
- filter_function(function_name: str) Parameters[source]
Filters based on the name of the
InteractionFunctionof eachParameter- Parameters:
function_name (str) – The name of the
InteractionFunctionofParameterobjects to return, for example'LennardJones'or'HarmonicPotential'.- Returns:
The
Parameterobjects which have afunctionwith the specifiedfunction_name- Return type:
- filter_interaction(interaction_name: str) Parameters[source]
Filters based on the name of the
Interactionof eachParameter- Parameters:
interaction_name (str) – The name of the
InteractionofParameterobjects to return, for example'Bond'.- Returns:
The
Parameterobjects which have anInteractionwith the specifiedinteraction_name- Return type:
- filter_name(name: str) Parameters[source]
Filters by
name- Parameters:
name (str) – The
nameof theParameterobjects to return.- Returns:
The
Parameterobjects withname- Return type:
- filter_structure(structure_name: str) Parameters[source]
Filters based on the name of the
Structureto which eachParameterapplies- Parameters:
structure_name (str) – The name of a
Structure.- Returns:
The
Parameterobjects which are applied to aStructurewhich has the specifiedstructure_name- Return type:
- filter_value(comparison: str, value: float) Parameters[source]
Filters by
value- Parameters:
- Returns:
The
Parameterobjects which return a True when their values are compared withvalueusing thecomparisonoperator- Return type:
- class MDMC.MD.Periodic(K1: float, n1: float, d1: float, *parameters: float)[source]
Bases:
InteractionFunctionPeriodic 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^-1n1 (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
Periodicfor aDihedralAngleinteraction withimproper==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
Periodicfor aDihedralAngleinteraction withimproper==False(i.e. a proper dihedral), withK1=3.53548,K2=-4.02501andK3=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)
- class MDMC.MD.Rattle(accuracy: float, max_iterations: int)[source]
Bases:
ConstraintAlgorithmHolds the parameters which are required for the RATTLE algorithm to be applied to the constrained interactions
- class MDMC.MD.Shake(accuracy: float, max_iterations: int)[source]
Bases:
ConstraintAlgorithmHolds the parameters which are required for the SHAKE algorithm to be applied to the constrained interactions
- class MDMC.MD.Simulation(universe: Universe, traj_step: int, time_step: float = 1.0, engine: str = 'lammps', **settings)[source]
Bases:
objectSets 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
Universeon which the simulation is performed.traj_step (int) – How many steps the simulation should take between dumping each
CompactTrajectoryframe. Along withtime_stepdetermines the time separation of calculated variables such as energy.time_step (float, optional) – Simulation timestep in
fs. Default is 1.engine (str, optional) – The
MDEngineused 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 thetemperature. If all atoms in theuniversehave 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 theengine, not theuniverse, are affected.integrator(str)Simulation time integrator.
lj_options(dict)option:valuepairs for Lennard-Jones interactions.es_options(dict)option:valuepairs 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.
- traj_step
How many steps the simulation should take between dumping each
CompactTrajectoryframe. Along withtime_stepdetermines the time separation of calculated variables such as energy.- Type:
- engine
A subclass of
MDEnginewhich provides the interface to the MD library. Default is'lammps'.- Type:
MDEngine, optional
- time_step
- trajectory
- auto_equilibrate(variables: list[str] = ('temp', 'pe'), eq_step: int = 10, window_size: int = 100, tolerance: float = 0.01, max_step=10000) 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
equilibrationis 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
CompactTrajectorybeing stored- Return type:
int
- property trajectory: CompactTrajectory | None
The
CompactTrajectoryproduced by the most recent production run of theSimulation.- Returns:
Most recent production run
CompactTrajectory, or None if no production run has been performed- Return type:
- class MDMC.MD.Structure(position: list[float] | tuple[float] | ndarray, velocity: list[float] | tuple[float] | ndarray, name: str)[source]
Bases:
ABCAbstract base class for all structural units
- Parameters:
position (list, tuple, numpy.ndarray) – A 3 element list, tuple or
arraywith the position in units ofAng.velocity (list, tuple, numpy.ndarray) – A 3 element list, tuple or
arraywith the velocity in units ofAng / fs.name (str) –
The name of the structure.
Attributes
---------- – ID : int A unique identifier for each
Structure.universe (Universe) – The
Universeto which theStructurebelongs.name – The name of the structure.
parent (Structure) –
Structureto which this unit belongs, orself
- property atoms: list[Atom]
Get a list of all of the Atom objects in the structure by recursively calling
atomsfor all substructures- Returns:
All atoms in the structure
- Return type:
- abstract property bonded_interaction_pairs: list[tuple[Interaction, tuple[Atom]]]
Get bonded interactions acting on the
Structureand the other atoms to which the atom is bonded- Returns:
(
interaction,atoms) pairs acting on theStructure, whereatomsis a tuple of allAtomobjects for that specific bondedinteraction. At least one of theseAtomobjects belongs to theStructure- Return type:
Example
For an O
Atomwith 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:
BondedInteractionsacting on theStructure- Return type:
- 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
arrayof float specifying thepositionof the newStructure- Returns:
A
Structureof the same type with all non-unique attributes copied and a newposition- Return type:
- property interactions: list[Interaction]
Get a list of the interactions acting on the
Structure- Returns:
Interactions acting on the
Structure- Return type:
- abstract is_equivalent(structure: Structure) bool[source]
Checks the passed
Structureagainst 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:
- abstract property nonbonded_interactions: list[NonBondedInteraction]
Get a list of the nonbonded interactions acting on the
Structure- Returns:
NonBondedInteractionsacting on theStructure- Return type:
- property position: ndarray
Get or set the position of the center of mass of the
StructureinAng- Return type:
- property structure_type: str
Get the class of the
Structure.- Returns:
The name of the class
- Return type:
- property top_level_structure: Structure
Get the top level structure (i.e.
Structurewhich has noparent) of theStructure- Returns:
Highest level
Structureof which it is a member- Return type:
- translate(displacement: tuple | ndarray) None[source]
Translate the structural unit by the specified displacement
- Parameters:
Displacement (tuple, numpy.ndarray) – A three element tuple or
arrayof float
- abstract property universe: 'Universe' | None
Get the
Universeto which theStructurebelongs- Returns:
The
Universeto which theStructurebelongs or None- Return type:
- valid_position(position: list[float] | tuple[float] | ndarray = None) bool[source]
Checks if the specified
positionis within the bounds of theStructure.universe, if it has one- Parameters:
position (list, tuple, numpy.ndarray) – 3 element list, tuple or
arraywith units ofAngor None. If None then thepositionof theStructureis used.- Returns:
True if
positionis withinUniverseor there is no associatedUniverse. False ifStructurehas an associatedUniversebut thepositionis not within its bounds.- Return type:
- Raises:
ValueError – If
positionif undefined
- class MDMC.MD.Universe(dimensions, force_field=None, structures=None, **settings)[source]
Bases:
AtomContainerClass where the configuration and topology are defined
- Parameters:
dimensions (numpy.ndarray, list, float) – Dimensions of the
Universe, in units ofAng. 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: OPLSAA, TIP3P, TIP3PFB, SPCE, SPC. Default is None.
structures (list, optional) –
Structureobjects contained in theUniverse. 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_solverordispersive_solvermay 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
Universein units ofAng.- Type:
- configuration
Stores the content, i.e. configuration of atoms etc within the universe
- Type:
- 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:
- electrostatic_solver
The k-space solver to be used for electrostatic interactions.
- Type:
- dispersive_solver
The k-space solver to be used for dispersive interactions.
- Type:
- constraint_algorithm
The constraint algorithm which will be applied to constrained
BondedInteractions.- Type:
- 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, whereatomsis a tuple of all atoms for that specific bondedinteraction
- add_force_field(force_field: str, *interactions: Interaction, **settings: dict) None[source]
Adds a force field to the specified
interactions. If nointeractionsare passed, the force field is applied to all interactions in theUniverse.- Parameters:
force_field (str) – The
ForceFieldto parameterizeinteractions(if provided), or all theinteractionsin theUniverse. The availableForceFieldare: OPLSAA, TIP3P, TIP3PFB, SPCE, SPC*interactions –
Interactionobjects to parameterize with theForceField**settings –
add_dispersions(bool or list ofAtoms)If True, a
Dispersioninteraction will be added to all atoms in theUniverse. If a list ofAtomobjects is provided, theDispersionwill be added to these instead. Any addedDispersioninteractions (and any previously defined) will then be parametrized by theForceField. TheDispersioninteractions added will only be like-like. By default, noDispersioninteractions 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
Structureto theUniverse, with optionalForceFieldapplying only to thatStructure- Parameters:
structure (Structure or int) – The
Structure(or itsIDas an int) to be added to theUniverseforce_field (str, optional) – The force field to be applied to the structure. The available
ForceFieldare: OPLSAA, TIP3P, TIP3PFB, SPCE, SPCcenter (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:interactionspairs whereatom_typeis a int specifying the atom type andinteractionsis a list ofInteractionobjects acting on thatatom_type- Return type:
- property atom_types: list[Atom]
Get the atom types of atoms in the
Universe- Returns:
The atom types in the
Universe- Return type:
- property atoms: list[Atom]
Get a list of the atoms in the
Universe- Returns:
The atoms in the
Universe- Return type:
- property bonded_interaction_pairs: list
Get the bonded interactions and the atoms they apply to
- Returns:
The (
interaction,atoms) pairs in theUniverse, whereatomsis a tuple of allAtomobjects for that specificinteraction- Return type:
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:
- property density
- property dimensions
- property element_dict: dict[str, Atom]
Get the elements in the
Universeand exampleAtomobjects for each elementThis is required for MD engines which assign the same potential parameters for all identical element types.
- Returns:
element:atom pairs, whereatomis a singleAtomof the specifiedelement.- Return type:
- property element_list: list[str]
The elements of the atoms in the
Universe- Returns:
The elements in the
Universe- Return type:
- property element_lookup: dict[str, Atom]
Get the elements by ID in the
UniverseThis is required for MD engines which assign the same potential parameters for all identical element names.
- Returns:
element:atom pairs, whereatomis a singleAtomof the specifiedelement.- Return type:
- property equivalent_top_level_structures_dict: dict[Structure, int]
Get a dict of equivalent top level
Structureobjects that exist in theUniverseas keys, and the number ofStructurethat are equivalent to the key as a value. This does not include anyStructurethat are a subunit of another structure belonging to theUniverse.- Returns:
The top level
Structureobjects in theUniverseas keys, and the number of each as a value- Return type:
- fill(structures: Structure, force_field: str = None, num_density: float = None, num_struc_units: int = None) None[source]
A liquid-like filling of the
Universeindependent of existing atomsAdds copies of
structuresto existing configuration untilUniverseis 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
Structurewith which to fill theUniverseforce_field (str) – Applies a
ForceFieldto theUniverse. The availableForceFieldare: OPLSAA, TIP3P, TIP3PFB, SPCE, SPCnum_density (float) – Non-negative float specifying the number density of the
Structure, in units ofStructure / Ang ^ -3num_struc_units (int) – Non-negative int specifying the number of passed
Structureobjects that the universe should be filled with, regardless ofUniverse.dimensions.
- Raises:
ValueError – If both
num_densityandnum_struc_unitsare passedValueError – If neither
num_densityornum_struc_unitsare passed
- property force_fields: None
Get the
ForceFieldacting on theUniverseThe available force fields are: OPLSAA, TIP3P, TIP3PFB, SPCE, SPC
- Returns:
ForceFieldthat applies to theUniverse- Return type:
- property interactions: list
Get the interactions in the
Universe- Returns:
The interactions in the
Universe- Return type:
- property molecule_list: list[Molecule]
Get a list of the
Moleculeobjects in theUniverse- Returns:
The
Moleculeobjects in theUniverse- Return type:
- property n_atoms: int
Get the number of atoms in the
Universe- Returns:
The number of atoms in the
Universe- Return type:
- property n_bonded: int
Get the number of bonded interactions in the
Universe- Returns:
The number of bonded interactions in the
Universe- Return type:
- property n_interactions: int
Get the number of interactions in the
Universe- Returns:
The number of interactions in the
Universe- Return type:
- property n_molecules: int
Get the number of molecules in the
Universe- Returns:
The number of molecules in the
Universe- Return type:
- property n_nonbonded: int
Get the number of nonbonded interactions in the
Universe- Returns:
The number of nonbonded interactions in the
Universe- Return type:
- property nbis_by_atom_type_pairs: dict[tuple[int], list[Interaction]]
Generates a dict of all nonbonded interactions possessed by all combinations of
atom_typepairs in theUniverse.- Returns:
- A dict of
{pair: interactions}pairs, where: pairis a tuple for a pair ofatom_typesin theUniverseinteractionsis a list of nonbonded interactions that exist for thispairofatom_types
Any
Dispersionsininteractionsare ones that exist explicity for thispair, whereas anyCoulombicsininteractionsare ones that exist for either of theatom_typesinpair.- A dict of
- Return type:
- property nonbonded_interactions: list
Get the nonbonded interactions in the
Universe- Returns:
The nonbonded interactions in the
Universe- Return type:
- property parameters: Parameters
Get the parameters of the interactions that exist within the
Universe- Returns:
The
Parametersobjects defined withinUniverse- Return type:
- solvate(density: float, tolerance: float = 1.0, solvent: str = 'SPCE', max_iterations: int = 100, **settings: dict) None[source]
Fills the
Universewith solvent molecules according to pre-defined coordinates.- Parameters:
density (float) – The desired density of the
Solventthat solvates theUniverse, in units ofamu Ang ^ -3tolerance (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
Solventfrom 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
ConstraintAlgorithmwhich is applied to theUniverse. If an inbuiltSolventis selected (e.g. ‘SPCE’) andconstraint_algorithmis not passed, theConstraintAlgorithmwill default toShake(1e-4, 100).
- Raises:
ValueError – If the
Universehas already been solvated with a different density.
- property solvent_density
- property structure_list: None
Get a list of all
Structureobjects that exist in theUniverse. This includes allStructurethat are a subunit of another structure belonging to theUniverse.- Returns:
The
Structureobjects in theUniverse- Return type:
- property top_level_structure_list: list[Structure]
Get a list of the top level
Structureobjects that exist in theUniverse. This does not include anyStructurethat are a subunit of another structure belonging to theUniverse.- Returns:
The top level
Structureobjects in theUniverse- Return type:
- property volume
- MDMC.MD.filter_atoms(atoms: list[Atom], predicate: Callable) list[Atom][source]
Filters a list of Atoms with a given predicate
- MDMC.MD.filter_atoms_element(atoms: list[Atom], element: str) list[Atom][source]
Filters a list of atoms based on the atomic element
- MDMC.MD.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 onlyn_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 ofsymbols.- Return type:
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'
- MDMC.MD.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 ofInteractionFunction.- 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 (orUnit) 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 inparameter_units.
Examples
The following adds units of
'Ang'to parameteralpha, units of's'to the parameterbeta, and units of'atm'to the parametergamma:@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
UnitFloatorUnitNDArraywill not be created). So to setepsilonas unitless:@inter_func_decorator(('delta', 'arb'), ('epsilon', None), ('gamma', 'deg')) def __init__(self, delta, epsilon, gamma): ...