pub struct MdState {
pub cfg: MdConfig,
pub atoms: Vec<AtomDynamics>,
pub adjacency_list: Vec<Vec<usize>>,
pub force_field_params: ForceFieldParamsIndexed,
pub time: f64,
pub step_count: usize,
pub snapshots: Vec<Snapshot>,
pub cell: SimBox,
pub neighbors_nb: NeighborsNb,
pub water_pme_sites_forces: Vec<[Vec3; 3]>,
pub kinetic_energy: f64,
pub potential_energy: f64,
/* private fields */
}
Fields§
§cfg: MdConfig
§atoms: Vec<AtomDynamics>
§adjacency_list: Vec<Vec<usize>>
§force_field_params: ForceFieldParamsIndexed
§time: f64
Current simulation time, in picoseconds.
step_count: usize
§snapshots: Vec<Snapshot>
These are the snapshots we keep in memory, accumulating.
cell: SimBox
§neighbors_nb: NeighborsNb
§water_pme_sites_forces: Vec<[Vec3; 3]>
§kinetic_energy: f64
kcal/mol
potential_energy: f64
Implementations§
Source§impl MdState
impl MdState
pub fn apply_bond_stretching_forces(&mut self)
Sourcepub fn apply_angle_bending_forces(&mut self)
pub fn apply_angle_bending_forces(&mut self)
This maintains bond angles between sets of three atoms as they should be from hybridization. It reflects this hybridization, steric clashes, and partial double-bond character. This identifies deviations from the ideal angle, calculates restoring torque, and applies forces based on this to push the atoms back into their ideal positions in the molecule.
Valence angles, which are the angle formed by two adjacent bonds ba et bc in a same molecule; a valence angle tends to maintain constant the anglê abc. A valence angle is thus concerned by the positions of three atoms.
Sourcepub fn apply_dihedral_forces(&mut self, improper: bool)
pub fn apply_dihedral_forces(&mut self, improper: bool)
This maintains dihedral angles. (i.e. the angle between four atoms in a sequence). This models effects such as σ-bond overlap (e.g. staggered conformations), π-conjugation, which locks certain dihedrals near 0 or τ, and steric hindrance. (Bulky groups clashing).
This applies both “proper” linear dihedral angles, and “improper”, hub-and-spoke dihedrals. These two angles are calculated in the same way, but the covalent-bond arrangement of the 4 atoms differs.
Sourcepub fn shake_hydrogens(&mut self)
pub fn shake_hydrogens(&mut self)
Part of our SHAKE + RATTLE algorithms for fixed hydrogens.
Sourcepub fn rattle_hydrogens(&mut self, mol_i: usize)
pub fn rattle_hydrogens(&mut self, mol_i: usize)
Part of our SHAKE + RATTLE algorithms for fixed hydrogens.
Source§impl MdState
impl MdState
Sourcepub fn step(&mut self, dev: &ComputationDevice, dt: f32)
pub fn step(&mut self, dev: &ComputationDevice, dt: f32)
One Velocity-Verlet step (leap-frog style) of length dt
is in picoseconds (10^-12),
with typical values of 0.001, or 0.002ps (1 or 2fs).
This method orchestrates the dynamics at each time step.
Source§impl MdState
impl MdState
Sourcepub fn build_neighbors_if_needed(&mut self)
pub fn build_neighbors_if_needed(&mut self)
Call during each step; determines if we need to rebuild neighbors, and does so A/R. todo: Run on GPU?
Sourcepub fn rebuild_dy_water_inv(&mut self)
pub fn rebuild_dy_water_inv(&mut self)
This inverts our neighbor set between water and dynamic atoms.
Source§impl MdState
impl MdState
Sourcepub fn apply_nonbonded_forces(&mut self, dev: &ComputationDevice)
pub fn apply_nonbonded_forces(&mut self, dev: &ComputationDevice)
Applies Coulomb and Van der Waals (Lennard-Jones) forces on dynamic atoms, in place. We use the MD-standard [S]PME approach to handle approximated Coulomb forces. This function applies forces from dynamic, static, and water sources.
Source§impl MdState
impl MdState
Sourcepub fn water_vv_first_half_and_drift(&mut self, dt: f32, dt_half: f32)
pub fn water_vv_first_half_and_drift(&mut self, dt: f32, dt_half: f32)
Verlet velocity integration for water, part 1. Forces for this step must be pre-calculated. Accepts as mutable to allow projecting M/EP force onto the other atoms.
In addition to the VV half-kick and drift, it handles force projection from M/EP, and applying SETTLE to main each molecul’s rigid geometry.
Sourcepub fn water_vv_second_half(&mut self, dt_half: f32)
pub fn water_vv_second_half(&mut self, dt_half: f32)
Velocity-Verlet integration for water, part 2. Forces (as .accel) must be computed prior to this step.
Source§impl MdState
impl MdState
pub fn new( cfg: &MdConfig, mols: &[MolDynamics<'_>], param_set: &FfParamSet, ) -> Result<Self, ParamError>
pub fn save_snapshots_to_file(&mut self, path: &Path)
Trait Implementations§
Auto Trait Implementations§
impl Freeze for MdState
impl !RefUnwindSafe for MdState
impl !Send for MdState
impl !Sync for MdState
impl Unpin for MdState
impl !UnwindSafe for MdState
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left
is true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left(&self)
returns true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§impl<T> Pointable for T
impl<T> Pointable for T
Source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
Source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self
from the equivalent element of its
superset. Read moreSource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self
is actually part of its subset T
(and can be converted to it).Source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset
but without any property checks. Always succeeds.Source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self
to the equivalent element of its superset.