[]Struct lumol::sys::System

pub struct System {
    pub simulated_degrees_of_freedom: DegreesOfFreedom,
    pub step: u64,
    // some fields omitted
}

The System type hold all the data about a simulated system.

This data contains:

In the implementation, the particles contained in a molecule are guaranteed to be contiguous in memory. This allow for faster access when iterating over molecules, and easier molecule removal from the system.

Fields

simulated_degrees_of_freedom: DegreesOfFreedom

Number of degrees of freedom simulated in the system. This default to DegreesOfFreedom::Particles, and is set in the simulation setup.

step: u64

The current simulation step

Implementations

impl System

pub fn new() -> System

Create a new empty System

pub fn with_cell(cell: UnitCell) -> System

Create an empty system with a specific unit cell

pub fn add_molecule(&mut self, molecule: Molecule)

Add a molecule to the system

pub fn composition(&self) -> Composition

Get the composition in particles and molecules of the configuration

pub fn simulated_temperature(&mut self, temperature: Option<f64>)

Use an external temperature for all the system properties. Calling this with Some(temperature) will replace all the computation of the temperature from the velocities with the given values. Calling it with None will use the velocities.

The default is to use the velocities unless this function is called.

impl System

Functions related to interactions

pub fn energy_evaluator(&self) -> EnergyEvaluator<'_>

Get an helper struct to evaluate the energy of this system.

pub fn set_pair_potential(&mut self, (&str, &str), potential: PairInteraction)

Set the pair interaction potential for atoms with types i and j

pub fn set_bond_potential(
    &mut self,
    (&str, &str),
    potential: Box<dyn BondPotential + 'static, Global>
)

Set the bond interaction potential for atoms with types i and j

pub fn set_angle_potential(
    &mut self,
    (&str, &str, &str),
    potential: Box<dyn AnglePotential + 'static, Global>
)

Set the angle interaction potential for atoms with types i, j, and k

pub fn set_dihedral_potential(
    &mut self,
    (&str, &str, &str, &str),
    potential: Box<dyn DihedralPotential + 'static, Global>
)

Set the dihedral angle interaction potential for atoms with types i, j, k, and m.

pub fn set_coulomb_potential(
    &mut self,
    potential: Box<dyn CoulombicPotential + 'static, Global>
)

Set the coulombic interaction for all pairs to potential

pub fn add_global_potential(
    &mut self,
    potential: Box<dyn GlobalPotential + 'static, Global>
)

Add the potential global interaction

pub fn pair_potential(&self, i: usize, j: usize) -> Option<&PairInteraction>

Get the pair potential acting between the particles at indexes i and j.

pub fn bond_potential(&self, i: usize, j: usize) -> Option<&dyn BondPotential>

Get the bond potential acting between the particles at indexes i and j.

pub fn angle_potential(
    &self,
    i: usize,
    j: usize,
    k: usize
) -> Option<&dyn AnglePotential>

Get the angle potential acting between the particles at indexes i, j and k.

pub fn dihedral_potential(
    &self,
    i: usize,
    j: usize,
    k: usize,
    m: usize
) -> Option<&dyn DihedralPotential>

Get the dihedral angles potential acting between the particles at indexes i, j, k and m.

pub fn coulomb_potential(&self) -> Option<&dyn CoulombicPotential>

Get the coulombic interaction for the system

pub fn global_potentials(&self) -> &[Box<dyn GlobalPotential + 'static, Global>]

Get all global interactions for the system

pub fn maximum_cutoff(&self) -> Option<f64>

Get maximum cutoff from coulomb, pairs and global interactions.

impl System

Functions to get physical properties of a system.

pub fn degrees_of_freedom(&self) -> usize

Get the number of degrees of freedom in the system

pub fn kinetic_energy(&self) -> f64

Get the kinetic energy of the system.

pub fn potential_energy(&self) -> f64

Get the potential energy of the system.

pub fn total_energy(&self) -> f64

Get the total energy of the system.

pub fn temperature(&self) -> f64

Get the temperature of the system.

pub fn volume(&self) -> f64

Get the volume of the system.

pub fn virial(&self) -> Matrix3

Get the virial of the system as a tensor

pub fn pressure(&self) -> f64

Get the pressure of the system from the virial equation, at the system instantaneous temperature.

pub fn stress(&self) -> Matrix3

Get the stress tensor of the system from the virial equation.

pub fn forces(&self) -> Vec<Vector3D>

Get the forces acting on all the particles in the system

impl System

pub fn check(&self)

Check the system before running a simulation

Methods from Deref<Target = Configuration>

pub fn are_in_same_molecule(&self, i: usize, j: usize) -> bool

Check if the particles at indexes i and j are in the same molecule

pub fn molecules(&self) -> MoleculeIter<'_>

Notable traits for MoleculeIter<'a>

impl<'a> Iterator for MoleculeIter<'a> type Item = MoleculeRef<'a>;

Get an iterator over the molecules in the configuration.

pub fn molecules_mut(&mut self) -> MoleculeIterMut<'_>

Notable traits for MoleculeIterMut<'a>

impl<'a> Iterator for MoleculeIterMut<'a> type Item = MoleculeRefMut<'a>;

Get an iterator over the molecules in the configuration.

pub fn molecule(&self, id: usize) -> MoleculeRef<'_>

Get the molecule at index id

pub fn molecule_mut(&mut self, id: usize) -> MoleculeRefMut<'_>

Get the molecule at index id

pub fn molecule_id(&self, i: usize) -> usize

Get the index of the molecule containing the particle i

pub fn bond_path(&self, i: usize, j: usize) -> BondPath

Get the length of the shortest bond path to go from the particle i to the particle j. If the particles are not in the same molecule, the length is -1. Else, this length is 0 if i == j, 1 if there is a bond between i and j, etc.

pub fn remove_molecule(&mut self, molid: usize)

Remove the molecule at index i

pub fn add_bond(
    &mut self,
    particle_i: usize,
    particle_j: usize
) -> Vec<Permutation>

Add a bond between the particles at indexes i and j. The particles should have been added to the configuration before calling this.

Warning

If the bond is between two particles which are not in the same molecule, the two molecules are merged together by moving particles in the particles list, and thus invalidate any previously stored index. In particular, any bond, angle, dihedral or molecule is invalidated.

This function will return the list of atomic permutations that where applied in order to ensure that molecules are contiguous in memory.

pub fn add_molecule(&mut self, molecule: Molecule)

Add a molecule to the configuration, putting the new particles at the end of the particles list

pub fn size(&self) -> usize

Get the number of particles in this configuration

pub fn is_empty(&self) -> bool

Check if this configuration contains any particle

pub fn center_of_mass(&self) -> Vector3D

Return the center-of-mass of the configuration

pub fn particles(&self) -> ParticleSlice<'_>

Get the list of particles in this configuration, as a ParticleSlice.

pub fn particles_mut(&mut self) -> ParticleSliceMut<'_>

Get the list of particles in this configuration, as a mutable ParticleSliceMut.

pub fn distance(&self, i: usize, j: usize) -> f64

Get the distance between the particles at indexes i and j

pub fn nearest_image(&self, i: usize, j: usize) -> Vector3D

Get the vector between the nearest image of particle j with respect to particle i.

pub fn angle(&self, i: usize, j: usize, k: usize) -> f64

Get the angle between the particles i, j and k

pub fn angle_and_derivatives(
    &self,
    i: usize,
    j: usize,
    k: usize
) -> (f64, Vector3D, Vector3D, Vector3D)

Get the angle and the derivatives of the angle between the particles i, j and k

pub fn dihedral(&self, i: usize, j: usize, k: usize, m: usize) -> f64

Get the dihedral angle between the particles i, j, k and m

pub fn dihedral_and_derivatives(
    &self,
    i: usize,
    j: usize,
    k: usize,
    m: usize
) -> (f64, Vector3D, Vector3D, Vector3D, Vector3D)

Get the dihedral angle and the derivatives of the dihedral angle between the particles i, j, k and m

Trait Implementations

impl Clone for System

impl Deref for System

type Target = Configuration

The resulting type after dereferencing.

impl DerefMut for System

impl From<Frame> for System

Auto Trait Implementations

impl !RefUnwindSafe for System

impl Send for System

impl Sync for System

impl Unpin for System

impl !UnwindSafe for System

Blanket Implementations

impl<T> Any for T where
    T: 'static + ?Sized
[src]

impl<T> Borrow<T> for T where
    T: ?Sized
[src]

impl<T> BorrowMut<T> for T where
    T: ?Sized
[src]

impl<T> From<T> for T[src]

impl<T, U> Into<U> for T where
    U: From<T>, 
[src]

impl<T> Pointable for T

type Init = T

The type for initializers.

impl<T> ToOwned for T where
    T: Clone
[src]

type Owned = T

The resulting type after obtaining ownership.

impl<T, U> TryFrom<U> for T where
    U: Into<T>, 
[src]

type Error = Infallible

The type returned in the event of a conversion error.

impl<T, U> TryInto<U> for T where
    U: TryFrom<T>, 
[src]

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.

impl<V, T> VZip<V> for T where
    V: MultiLane<T>,