MDMC.trajectory_analysis package
Subpackages
- MDMC.trajectory_analysis.observables package
- Submodules
- MDMC.trajectory_analysis.observables.concurrency_tools module
- MDMC.trajectory_analysis.observables.fqt module
AbstractFQt
AbstractFQt.FQt
AbstractFQt.apply_resolution()
AbstractFQt.calculate_SQw()
AbstractFQt.calculate_from_MD()
AbstractFQt.calculate_rho_config()
AbstractFQt.dependent_variables
AbstractFQt.dependent_variables_structure
AbstractFQt.errors
AbstractFQt.independent_variables
AbstractFQt.t
AbstractFQt.uniformity_requirements
FQt
calc_incoherent_scatt_length()
calculate_rho()
get_point_group()
wyckoff_symmetries()
- MDMC.trajectory_analysis.observables.fqt_coh module
- MDMC.trajectory_analysis.observables.fqt_incoh module
- MDMC.trajectory_analysis.observables.obs module
Observable
Observable.reader
Observable.calculate_from_MD()
Observable.data
Observable.dependent_variables
Observable.dependent_variables_structure
Observable.errors
Observable.independent_variables
Observable.maximum_frames()
Observable.minimum_frames()
Observable.name
Observable.origin
Observable.read_from_file()
Observable.uniformity_requirements
Observable.use_FFT
- MDMC.trajectory_analysis.observables.obs_factory module
- MDMC.trajectory_analysis.observables.pdf module
PairDistributionFunction
PairDistributionFunction.PDF
PairDistributionFunction.PDF_err
PairDistributionFunction.calculate_from_MD()
PairDistributionFunction.dependent_variables
PairDistributionFunction.dependent_variables_structure
PairDistributionFunction.errors
PairDistributionFunction.independent_variables
PairDistributionFunction.maximum_frames()
PairDistributionFunction.minimum_frames()
PairDistributionFunction.name
PairDistributionFunction.r
PairDistributionFunction.uniformity_requirements
- MDMC.trajectory_analysis.observables.sqw module
AbstractSQw
AbstractSQw.E
AbstractSQw.SQw
AbstractSQw.SQw_err
AbstractSQw.calculate_E()
AbstractSQw.calculate_dt()
AbstractSQw.calculate_from_MD()
AbstractSQw.calculate_resolution_functions()
AbstractSQw.dependent_variables
AbstractSQw.dependent_variables_structure
AbstractSQw.errors
AbstractSQw.independent_variables
AbstractSQw.recreated_Q
AbstractSQw.uniformity_requirements
AbstractSQw.validate_energy()
AbstractSQw.w
SQw
SQwCoherent
SQwIncoherent
SQwMixins
- Module contents
- Observables
CoherentDynamicStructureFactor
CoherentIntermediateScatteringFunction
DynamicStructureFactor
FQt
FQtCoh
FQtCoherent
FQtIncoherentFQtIncoh
FQt_coh
FQt_incoh
IncoherentDynamicStructureFactor
IncoherentIntermediateScatteringFunction
IntermediateScatteringFunction
PDF
PairDistributionFunction
PairDistributionFunction.PDF
PairDistributionFunction.PDF_err
PairDistributionFunction.calculate_from_MD()
PairDistributionFunction.dependent_variables
PairDistributionFunction.dependent_variables_structure
PairDistributionFunction.errors
PairDistributionFunction.independent_variables
PairDistributionFunction.maximum_frames()
PairDistributionFunction.minimum_frames()
PairDistributionFunction.name
PairDistributionFunction.r
PairDistributionFunction.uniformity_requirements
SQw
SQwCoh
SQwCoherent
SQwIncoherentSQwIncoh
SQw_coh
SQw_incoh
Submodules
MDMC.trajectory_analysis.compact_trajectory module
Module for memory-efficient MD trajectory handling.
Notes
Seeing how the current (as of August 2022) trajectory implementation consumes a lot of memory, an attempt is being made to build an object that will contain the bare minimum of functionality, to show the limits of performance that we can achieve within Python.
- class MDMC.trajectory_analysis.compact_trajectory.CompactTrajectory(n_steps: int = 0, n_atoms: int = 1, useVelocity: bool = False, bits_per_number: int = 8, **settings: dict)[source]
Bases:
object
Store an MD trajectory in numpy arrays.
The units are system units.
The normal workflow for a CompactTrajectory is:
Create and allocate memory.
traj = CompactTrajectory(n_step, n_atoms, useVelocity) or traj = CompactTrajectory() traj.preAllocate(n_step, n_atoms, useVelocity)
Write values into arrays.
for n in frame_numbers: if not traj.validateTypes(atom_types): traj.WriteOneStep(step_number = n, time = time_value, positions = pos_values) traj.setDimensions(box_dimensions, step_number = n)
Add chemical element information.
traj.labelAtoms(atom_symbols, atom_masses) traj.setCharges(atom_charges)
Check internal consistency, trim arrays.
traj.postProcess()
- Parameters:
n_steps (int, optional) –
Number of simulation steps in the trajectory. Defaults to 0.
n_steps = 0
means thatpreAllocate
will have to be run separately after theCompactTrajectory
has been created.n_steps > 0
means thatpreAllocate
will be run immediately in the__init__
function.n_atoms (int, optional) – Number of atoms in the system. Defaults to 1.
useVelocity (bool, optional) – If the trajectory contains velocities, set to
True
to allocate an additional array for the velocity values. Defaults to False.bits_per_number (int, optional) – Set number of bytes in floating point numbers.
**settings (dict) – Extra options.
- property data: ndarray
The step numbers and time steps in a single array.
- Returns:
Concatenation of separate data arrays.
- Return type:
- exportAtom(step_number: int = 0, atom_number: int = 0)[source]
Create an
Atom
object for a chosen time step of the simulation and atom number.- Parameters:
- Returns:
A single
Atom
object.- Return type:
Notes
For compatibility with
Trajectory
.
- filter_by_element(elements: list[str]) CompactTrajectory [source]
Filter subtrajectory by chemical symbol.
- filter_by_time(start: float, end: float = None) CompactTrajectory [source]
Filter to within time defined by
start
andend
.Create another CompactTrajectory, containing only the frames with the time values equal to or larger than the start parameter, and smaller than the end parameter.
If only one time value is given, the function will create a CompactTrajectory with a single time step equal to start, or raise an error if start is not an element of the time array.
- Parameters:
- Returns:
A copy with
times
in half open interval defined bystart
andend
.- Return type:
- Raises:
ValueError – No valid MD frames found.
- filter_by_type(types: list[int]) CompactTrajectory [source]
Filter subtrajectory by atom ID.
- Parameters:
- Returns:
A
CompactTrajectory
containing only the atoms of the specified type.- Return type:
- labelAtoms(atom_symbols: dict[int, str] = None, atom_masses: dict[int, float] = None) bool [source]
Create internal atom references.
A list of chemical elements of all the atoms in a trajectory frame,
a set of all the chemical elements present in a trajectory,
a list of atom masses of all the atoms in a trajectory frame.
It is up to the engine facade to construct the input dictionaries containing the correct information about the chemical elements and masses.
- property positions: ndarray
The position array.
Added for those parts of code that use the plural form instead of singular.
- Returns:
Time samples array.
- Return type:
- postProcess()[source]
Finalise trajectory analysis.
Notes
This function should be called after the all the trajectory steps have been read. It will discard the unnecessary rows of the arrays, in case we had allocated too many due to some rounding error.
- preAllocate(n_steps: int = 1, n_atoms: int = 1, useVelocity: bool = False)[source]
Allocate array space for data.
Creates empty arrays for storing atom positions, velocities, time values of the simulation frames, and the simulation box dimensions.
The numpy empty object do not initialise elements until they have been accessed by a read or write operation. This means that allocating too many steps is generally safe and harmless, as the unused steps will be removed when the postProcess method is called.
- Parameters:
n_steps (int, optional) – Number of simulation steps in the trajectory. Defaults to 1.
n_atoms (int, optional) – Number of atoms in the system. Defaults to 1.
useVelocity (bool, optional) – If the trajectory contains velocities, set to True to allocate an additional array for the velocity values. Defaults to False.
- setBytesPerNumber(bytes_per_number: int = 8)[source]
Change the number of bytes per number in the arrays.
This is used for storing atom positions, velocities, the frame timestamps and simulation box dimensions.
The best approach is to set the correct value before populating the arrays, but it is still possible to change the data type using this function when the
CompactTrajectory
already contains some values.- Parameters:
bytes_per_number (int, optional) – If 8, the arrays will use
np.float64
to store the positions, velocities and time at each step. Will always be rounded up to the nearest multiple of 2.
See also
CompactTrajectory._get_dtype
Type getter.
- setCharge(charge_list: list[float] = None) bool [source]
Set the values of partial charge for each atom.
It assumes that the charges are constant within the simulation.
- setDimensions(frame_dimensions: ndarray = None, step_num: int = -1)[source]
Write the simulation box dimensions into the object header.
Additionally, if the simulation box dimensions change with time, it will keep track of the new dimensions at each step in the self.changing_dimensions object.
- Parameters:
frame_dimensions (numpy.ndarray) – 3 float numbers defining the size of the simulation box along x, y and z.
step_num (int) – The number of the simulation frame at which the frame_dimensions array was read.
- subtrajectory(start: int = 0, stop: int = -1, step: int = 1, atom_filter: list[int] = None)[source]
Slice a
CompactTrajectory
.Returns another
CompactTrajectory
instance, which contains the same header information, and a subset of the original trajectory steps.The arrays in the original trajectory will be sliced following the pattern: new = old[start:stop:step].
Optionally,
atom_filter
can be specified to choose only specific atoms from the trajectory.It is recommended not to use
subtrajectory
directly in this case, but to use thefilter_by_element
orfilter_by_type
methods, which will create theelement_filter
themselves.- Parameters:
- Returns:
A trajectory containing the same header and the same or fewer steps than the original.
- Return type:
See also
filter_by_element
Get components by species.
filter_by_type
Get components by ID.
- property times: ndarray
The time array.
Added for those parts of code that use the plural form instead of singular.
- Returns:
Time samples array.
- Return type:
- validateTypes(raw_atom_types: array) bool [source]
Check and set the array of atom types.
If previously run, verify the array of atom types from the new frame is the same as the original array of atom types.
If the atom types have changed during the simulation, we cannot process the results using the
CompactTrajectory
object, and the validation will returnFalse
. It is up to the engine facade to decide what to do ifvalidateTypes
returnsFalse
.- Parameters:
raw_atom_types (numpy.array) –
An array of all the atom types, sorted by the atom ID.
If these are numbers, like in a LAMMPS simulation, then they will be stored directly.
If
raw_atom_types
are strings, an internal method called_create_types_from_elements
will create numbers that correspond to the strings, so that every atom type has its own number.- Returns:
Whether the atom types are unchanged.
- Return type:
- property velocities: ndarray
The velocity array.
Added for those parts of code that use the plural form instead of singular.
- Returns:
Velocities array.
- Return type:
- writeEmptyStep(step_num: int = -1, time: float = -1.0)[source]
Advance the iterators without writing any data.
It is basically a reduced version of the writeOneStep method. This was added since the tests/trajectory_analysis/test_PDF.py use a trajectory made of 1000 _empty_ configurations.
- writeOneStep(step_num: int = -1, time: float = -1.0, positions: array = None, velocities: array = None)[source]
Write the atom properties of a single frame into the trajectory arrays.
- Parameters:
step_num (int, optional) – The index at which the numbers will be written. Defaults to -1.
time (float) – The time stamp of the simulation step, in the correct time units (femtoseconds). Defaults to -1.0.
positions (np.array) – The array of the atom positions, shaped (n_atoms, 3). Defaults to None.
velocities (np.array, optional) – The array of the atom velocities, shaped (n_atoms, 3). If we don’t use velocities, it can be skipped. Defaults to None.
- MDMC.trajectory_analysis.compact_trajectory.configurations_as_compact_trajectory(*configs: list[TemporalConfiguration]) CompactTrajectory [source]
Populate
CompactTrajectory
arrays from list ofTemporalConfiguration
.- Parameters:
*configs (list[TemporalConfiguration]) – Most likely containing atoms.
- Returns:
With atom types and positions matching those in the input configurations.
- Return type:
- Raises:
TypeError – No configurations provided.
MDMC.trajectory_analysis.trajectory module
Module for Configuration
and related classes.
- class MDMC.trajectory_analysis.trajectory.AtomCollection[source]
Bases:
object
Base class for
Configurations
.
- class MDMC.trajectory_analysis.trajectory.Configuration(*structures: Structure, **settings: dict)[source]
Bases:
AtomCollection
A
Configuration
storesAtom
objects and their positions and velocities.- Parameters:
*structures – Zero or more
Structure
objects to be added to theConfiguration
.**settings (dict) – Extra options.
- Parameters:
universe (Universe) – The
Universe
of theConfiguration
.
- add_structure(structures: Structure) None [source]
Add the
Atom
objects from aStructure
to the data.- Parameters:
structures (Structure) – The
Structure
to add.
- property atom_positions: list[ndarray]
Get the list of
Atom.position
which belong to theConfiguration
.- Returns:
A list of
Atom.position
s inConfiguration
.- Return type:
- property atom_velocities: list[ndarray]
Get the list of
Atom.velocity
which belong to theConfiguration
.- Returns:
A list of
Atom.velocity` s in ``Configuration
.- Return type:
- property data: ndarray
Get or set the
Atom
properties of theConfiguration
.- Returns:
A structured NumPy
array
with'atom'
,'position'
, and'velocity'
fields.- Return type:
- property element_list: list[str]
Get the list of
Atom.element.symbol
which belong to theConfiguration
.
- element_set
- filter_atoms(predicate: Callable[[Atom], bool]) list[Atom] [source]
Filter the list of
Atom
using the predicate.
- filter_structures(predicate: Callable[[Structure], bool]) list[Structure] [source]
Filter the list of
Structures
using the predicate.
- scale(factor: float, vectors: Literal['positions', 'velocities'] = 'positions') None [source]
Scale either
atom_positions
oratom_velocities
by a factor.- Parameters:
factor (float) – Factor by which the vector is scaled.
vectors ({'positions', 'velocities'}) – Vectors to rescale.
- Raises:
NotImplementedError – THIS IS NOT IMPEMENTED.
Module contents
Modules relating to describing an MD trajectory and calculating observables.
Contents
observables
trajectory