physical_validation.data subpackage

physical_validation.data.simulation_data module

Data structures carrying simulation data.

class physical_validation.data.simulation_data.SimulationData(units=None, dt=None, system=None, ensemble=None, observables=None, trajectory=None)[source]

SimulationData: System information and simulation results

The SimulationData class holds both the information on the system and the results of a simulation run of that system. SimulationData contains all information on a simulation run needed by the physical validation tests. SimulationData objects can either be created directly by calling the class constructor, or by using a parser returning a SimulationData object.

static compatible(data_1, data_2) bool[source]

Checks whether two simulations are compatible for common validation.

Parameters
Returns

result

Return type

bool

property ensemble: physical_validation.data.ensemble_data.EnsembleData

Information on the sampled ensemble

Returns

ensemble

Return type

EnsembleData

Type

EnsembleData

property units: physical_validation.data.unit_data.UnitData

Information on the sampled units

Returns

units

Return type

UnitData

Type

UnitsData

property observables: physical_validation.data.observable_data.ObservableData

Observables collected during the simulation

Returns

observables

Return type

ObservableData

Type

ObservableData

property trajectory: physical_validation.data.trajectory_data.TrajectoryData

Trajectories collected during the simulation

Returns

trajectory

Return type

TrajectoryData

Type

TrajectoryData

property system: physical_validation.data.system_data.SystemData

Information on the system’s system

Returns

system

Return type

SystemData

Type

SystemData

property dt: float

The timestep of the simulation run.

Returns

timestep

Return type

float

set_ensemble(ensemble: str, natoms: Optional[float] = None, mu: Optional[float] = None, volume: Optional[float] = None, pressure: Optional[float] = None, energy: Optional[float] = None, temperature: Optional[float] = None) None[source]
raise_if_units_are_none(test_name: str, argument_name: str) None[source]

Raise if the unit data was not set

Parameters
  • test_name – String naming the test used for error output

  • argument_name – String naming the SimulationData argument used for error output

raise_if_ensemble_is_invalid(test_name: str, argument_name: str, check_pressure: bool, check_mu: bool) None[source]

Raise InputError if the ensemble data does not hold the required information needed for the tests.

Parameters
  • test_name – String naming the test used for error output

  • argument_name – String naming the SimulationData argument used for error output

  • check_pressure – Whether to check if the pressure is defined (NPT only).

  • check_mu – Whether to check if the chemical potential is defined (muVT only).

raise_if_system_data_is_invalid(test_name: str, argument_name: str, check_full_system_data_only: bool) None[source]

Raise InputError if the ensemble data does not hold the required information needed for the tests.

Parameters
  • test_name – String naming the test used for error output

  • argument_name – String naming the SimulationData argument used for error output

  • check_full_system_data_only – Whether to check the full system data only (number of atoms in the system, number of constraints in the system, number of translational and rotational degree of freedom reduction of the system). If false, this also checks the masses per atom and the molecule index and constraints.

raise_if_observable_data_is_invalid(required_observables: List[str], test_name: str, argument_name: str) None[source]

Raise InputError if there are missing observables entries, or if the required observables don’t have the same number of frames

Parameters
  • required_observables – List of strings denoting the observables required

  • test_name – String naming the test used for error output

  • argument_name – String naming the SimulationData argument used for error output

raise_if_trajectory_data_is_invalid(test_name: str, argument_name: str) None[source]

Raise InputError if there are missing observables entries, or if the required observables don’t have the same number of frames

Parameters
  • test_name – String naming the test used for error output

  • argument_name – String naming the SimulationData argument used for error output

physical_validation.data.unit_data module

Data structures carrying simulation data.

class physical_validation.data.unit_data.UnitData(kb: float, energy_conversion: float, length_conversion: float, volume_conversion: float, temperature_conversion: float, pressure_conversion: float, time_conversion: float, energy_str: str = 'ENE', length_str: str = 'LEN', volume_str: str = 'VOL', temperature_str: str = 'TEMP', pressure_str: str = 'PRESS', time_str: str = 'TIME')[source]

UnitData: Information about the units used

The information about units consists of different parts:

  • The name of the units (energy_str, length_str, volume_str,

    temperature_str, pressure_str, time_str),

  • the value of kB in the used energy units, and

  • the conversion factor to GROMACS units (kJ/mol, nm, nm^3, K, bar, ps).

The names are only used for output (console printing and plotting), and are optional. The conversion factors and kB are, on the other hand, used in computations and need to be given.

classmethod units(name: Optional[str] = None)[source]
property kb: float

The value of the Boltzmann constant

Type

float

property energy_str: str

Energy unit

Type

str

property length_str: str

Length unit

Type

str

property volume_str: str

Volume unit

Type

str

property temperature_str: str

Temperature unit

Type

str

property pressure_str: str

Pressure unit

Type

str

property time_str: str

Time unit

Type

str

property energy_conversion: float

Energy conversion factor, 1 energy_unit = energy_conversion * kJ/mol

Type

float

property length_conversion: float

Length conversion factor, 1 length_unit = length_conversion * nm

Type

float

property volume_conversion: float

Volume conversion factor, 1 volume_unit = volume_conversion * nm^3

Type

float

property temperature_conversion: float

Temperature conversion factor, 1 temperature_unit = temperature_conversion * K

Type

float

property pressure_conversion: float

Pressure conversion factor, 1 pressure_unit = pressure_conversion * bar

Type

float

property time_conversion: float

Time conversion factor, 1 time_unit = time_conversion * ps

Type

float

physical_validation.data.ensemble_data module

Data structures carrying simulation data.

class physical_validation.data.ensemble_data.EnsembleData(ensemble, natoms=None, mu=None, volume=None, pressure=None, energy=None, temperature=None)[source]

EnsembleData: Holds data defining the ensemble

The ensemble is a string indicating the thermodynamical ensemble a simulation was performed in, and is any of ‘NVE’, ‘NVT’, ‘NPT’, ‘muVT’. Depending on the ensemble, EnsembleData then holds additional information defining the ensemble, such as the number of particles N, the chemical potential mu, the volume V, the pressure P, the constant energy E or the temperature T. While any of these additional information are optional, most of them are needed by certain tests, such that not fully defining the ensemble results in warnings. The notable exception to this rule is the constant energy E for the NVE, which is not needed by any test and can hence be omitted without raising a warning.

static ensembles() Tuple[str, str, str, str][source]
property ensemble: str

Get ensemble

property natoms: int

Get natoms

property mu: numpy.ndarray

Get mu

property volume: float

Get volume

property pressure: float

Get pressure

property energy: float

Get energy

property temperature: float

Get temperature

physical_validation.data.trajectory_data module

Data structures carrying simulation data.

class physical_validation.data.trajectory_data.RectangularBox(box: numpy.ndarray)[source]
property box
gather(positions: numpy.ndarray, bonds: List[List[int]], molec_idx: List[int])[source]
class physical_validation.data.trajectory_data.TrajectoryData(position: Optional[Any] = None, velocity: Optional[Any] = None)[source]

TrajectoryData: The position and velocity trajectory along the simulation

The full trajectory is needed to calculate the equipartition of the kinetic energy. As they are used in connection, the position and velocity trajectories are expected to have the same shape and number of frames.

The position and velocity trajectories can be accessed either using the getters of an object, as in

  • trajectory.position

  • trajectory.velocity

or using the key notation, as in

  • trajectory[‘position’]

  • trajectory[‘velocity’]

static trajectories() Tuple[str, str][source]
property position: Optional[numpy.ndarray]

Get position

property velocity: Optional[numpy.ndarray]

Get velocity

physical_validation.data.observable_data module

Data structures carrying simulation data.

class physical_validation.data.observable_data.ObservableData(kinetic_energy: Optional[Any] = None, potential_energy: Optional[Any] = None, total_energy: Optional[Any] = None, volume: Optional[Any] = None, pressure: Optional[Any] = None, temperature: Optional[Any] = None, constant_of_motion: Optional[Any] = None, number_of_species: Optional[Any] = None)[source]

ObservableData: The trajectory of (macroscopic) observables during the simulation

Stores a number of different observables:

  • kinetic_energy: the kinetic energy of the system,

  • potential_energy: the potential energy of the system,

  • total_energy: the total energy of the system,

  • volume: the volume of the system box,

  • pressure: the pressure of the system,

  • temperature: the temperature of the system,

  • constant_of_motion: a constant of motion of the trajectory.

The observable trajectories can be accessed either using the getters of an object, as in

observables.kinetic_energy

or using the key notation, as in

observables[‘kinetic_energy’]

static observables() List[str][source]
property kinetic_energy: Optional[numpy.ndarray]

Get kinetic_energy

property potential_energy: Optional[numpy.ndarray]

Get potential_energy

property total_energy: Optional[numpy.ndarray]

Get total_energy

property volume: Optional[numpy.ndarray]

Get volume

property pressure: Optional[numpy.ndarray]

Get pressure

property temperature: Optional[numpy.ndarray]

Get temperature

property constant_of_motion: Optional[numpy.ndarray]

Get constant_of_motion

property number_of_species: Optional[numpy.ndarray]

Get number_of_species

property nframes: Optional[int]

Get number of frames

property kinetic_energy_per_molecule: Optional[numpy.ndarray]

Get kinetic_energy per molecule - used internally

physical_validation.data.system_data module

Data structure carrying information on the simulated system.

class physical_validation.data.system_data.SystemData(natoms: Optional[int] = None, nconstraints: Optional[float] = None, ndof_reduction_tra: Optional[float] = None, ndof_reduction_rot: Optional[float] = None, mass: Optional[numpy.ndarray] = None, molecule_idx: Optional[numpy.ndarray] = None, nconstraints_per_molecule: Optional[numpy.ndarray] = None)[source]

SystemData: Information about the atoms and molecules in the system.

The information stored in SystemData objects describes the atoms and molecules in the system as far as the physical validation tests need it.

The system is described in terms of

  • natoms: the total number of atoms in the system

  • nconstraints: the total number of constraints in the system

  • ndof_reduction_tra: global reduction of translational degrees of freedom (e.g. due to constraining the center of mass of the system)

  • ndof_reduction_rot: global reduction of rotational degrees of freedom (e.g. due to constraining the center of mass of the system)

The atoms are described in terms of

  • mass: a list of the mass of every atom in the system

The molecules are described by

  • molecule_idx: a list with the indices first atoms of every molecule (this assumes that the atoms are sorted by molecule)

  • nconstraints_per_molecule: a list with the number of constraints in every molecule

Only used internally:

  • ndof_per_molecule: a list with the number of degrees of freedom of every molecule

Reserved for future use:

  • bonds

  • constrained_bonds

Notes:

  • kinetic_energy.distribution() only requires information on the system (natoms, nconstraints, ndof_reduction_tra, ndof_reduction_rot)

  • kinetic_energy.equipartition() additionally requires information on the atoms and molecules (mass, molecule_idx, nconstraints_per_molecule)

All other tests do not require and information from SystemData.

property natoms: int

Number of atoms in the system

Type

int

property nconstraints: float

Total number of constraints in the system

Does not include the reduction of degrees of freedom in the absence of external forces.

Type

float

property ndof_reduction_tra: float

Number of translational degrees of freedom deducted from 3*[# of molecules]

Type

float

property ndof_reduction_rot: float

Number of rotational degrees of freedom deducted from 3*[# of molecules]

Type

float

property mass: numpy.ndarray

Mass vector for the atoms

Setter accepts array-like objects.

Type

nd-array

property molecule_idx: numpy.ndarray

List of index of first atom of each molecule

Setter accepts array-like objects.

Type

nd-array

property nconstraints_per_molecule: numpy.ndarray

List of number of constraints per molecule

Setter accepts array-like objects.

Type

nd-array

property ndof_per_molecule: Optional[List[Dict[str, float]]]

List of number of degrees of freedom per molecule

Type

nd-array

property bonds: List[List[int]]

List of bonds per molecule

Type

List[List[int]]

property constrained_bonds: List[List[int]]

List of constrained bonds per molecule

Type

List[List[int]]

physical_validation.data.parser module

Parsers read output files from MD simulation packages and create SimulationData objects with their contents.

class physical_validation.data.parser.Parser[source]

Parser base class

get_simulation_data() physical_validation.data.simulation_data.SimulationData[source]

physical_validation.data.gromacs_parser module

gromacs_parser.py

class physical_validation.data.gromacs_parser.GromacsParser(exe: Optional[str] = None, includepath: Optional[Union[str, List[str]]] = None)[source]
static units() physical_validation.data.unit_data.UnitData[source]
get_simulation_data(mdp: Optional[str] = None, top: Optional[str] = None, edr: Optional[str] = None, trr: Optional[str] = None, gro: Optional[str] = None) physical_validation.data.simulation_data.SimulationData[source]
Parameters
  • mdp (str, optional) – A string pointing to a .mdp file

  • top (str, optional) – A string pointing to a .top file

  • edr (str, optional) – A string pointing to a .edr file

  • trr (str, optional) – A string pointing to a .trr file

  • gro (str, optional) – A string pointing to a .gro file (Note: if also trr is given, gro is ignored)

Returns

result – A SimulationData filled with the results of the simulation as described by the provided GROMACS files.

Return type

SimulationData

physical_validation.data.flatfile_parser module

flatfile_parser.py

class physical_validation.data.flatfile_parser.FlatfileParser[source]
get_simulation_data(units: Optional[physical_validation.data.unit_data.UnitData] = None, ensemble: Optional[physical_validation.data.ensemble_data.EnsembleData] = None, system: Optional[physical_validation.data.system_data.SystemData] = None, dt: Optional[float] = None, position_file: Optional[str] = None, velocity_file: Optional[str] = None, kinetic_ene_file: Optional[str] = None, potential_ene_file: Optional[str] = None, total_ene_file: Optional[str] = None, volume_file: Optional[str] = None, pressure_file: Optional[str] = None, temperature_file: Optional[str] = None, const_of_mot_file: Optional[str] = None, number_of_species_file: Optional[str] = None) physical_validation.data.simulation_data.SimulationData[source]

Read simulation data from flat files

Returns a SimulationData object created from (optionally) provided UnitData, EnsembleData and SystemData, as well as TrajectoryData and ObservableData objects created from flat files. The files are expected to be in one of the following formats:

  • xyz-format trajectory files (position_file, velocity_file) - three numbers per line, separated by white space - frames delimited by a completely blank line - any character after (and including) a ‘#’ are ignored

  • 1d-format all other files - one number per line - any character after (and including) a ‘#’ are ignored

Parameters
  • units (UnitData, optional) – A UnitData object representing the units used in the simulation

  • ensemble (EnsembleData, optional) – A EnsembleData object representing the ensemble the simulation has been performed in

  • system (SystemData, optional) – A SystemData object representing the atoms and molecules in the system

  • dt (float, optional) – The time step used in the simulation

  • position_file (str, optional) – Path to a file in xyz-format containing the position trajectory

  • velocity_file (str, optional) – Path to a file in xyz-format containing the velocity trajectory

  • kinetic_ene_file (str, optional) – Path to a file in 1d-format containing the kinetic energy trajectory

  • potential_ene_file (str, optional) – Path to a file in 1d-format containing the potential energy trajectory

  • total_ene_file (str, optional) – Path to a file in 1d-format containing the total energy trajectory

  • volume_file (str, optional) – Path to a file in 1d-format containing the volume trajectory

  • pressure_file (str, optional) – Path to a file in 1d-format containing the pressure trajectory

  • temperature_file (str, optional) – Path to a file in 1d-format containing the temperature trajectory

  • const_of_mot_file (str, optional) – Path to a file in 1d-format containing the constant of motion trajectory

  • number_of_species_file (str, optional) – Path to a file in nd-format containing the number of species trajectory

Returns

result – A SimulationData filled with the provided ensemble and system objects as well as the trajectory data found in the edr and trr / gro files.

Return type

SimulationData