#![allow(non_upper_case_globals)]
#![allow(clippy::excessive_precision)]
use std::{
fmt,
fmt::{Display, Formatter},
};
#[allow(unused)]
#[cfg(target_arch = "x86_64")]
use lin_alg::f32::{Vec3x8, Vec3x16};
use lin_alg::{
f32::{Quaternion as QuaternionF32, Vec3 as Vec3F32, X_VEC, Z_VEC},
f64::Vec3,
};
use na_seq::Element;
use crate::{
AtomDynamics, KCAL_TO_NATIVE, MolDynamics, barostat::SimBox, non_bonded::CHARGE_UNIT_SCALER,
};
#[allow(unused)]
#[cfg(target_arch = "x86_64")]
use crate::{AtomDynamicsx8, AtomDynamicsx16};
pub(crate) mod init;
pub(crate) mod settle;
use settle::RA;
pub(crate) const O_MASS: f32 = 16.;
pub(crate) const H_MASS: f32 = 1.008;
pub(crate) const MASS_WATER_MOL: f32 = O_MASS + 2.0 * H_MASS;
pub(crate) const O_EP_R: f32 = 0.159_398_33;
pub(crate) const O_H_R: f32 = 0.872_433_13;
pub(crate) const H_O_H_θ: f32 = 1.808_161_105_066; const H_O_H_θ_HALF: f32 = 0.5 * H_O_H_θ;
const SIGMA_FACTOR: f32 = 2. / 1.122_462_048_309_373;
const O_RSTAR: f32 = 1.777_167_268;
pub const O_SIGMA: f32 = O_RSTAR * SIGMA_FACTOR;
pub const O_EPS: f32 = 0.212_800_813_0;
const Q_H: f32 = 0.6791 * CHARGE_UNIT_SCALER;
const Q_EP: f32 = -2. * Q_H;
const C_H: f32 = (O_EP_R / RA) / 2.;
const C_O: f32 = 1.0 - 2.0 * C_H;
pub(crate) const ACCEL_CONV_WATER_O: f32 = KCAL_TO_NATIVE / O_MASS;
pub(crate) const ACCEL_CONV_WATER_H: f32 = KCAL_TO_NATIVE / H_MASS;
#[derive(Clone, Debug, Default)]
pub enum Solvent {
None,
#[default]
WaterOpc,
WaterOpcSpecifyMolCount(usize),
Custom((Vec<(MolDynamics, usize)>, usize)),
}
impl Display for Solvent {
fn fmt(&self, f: &mut Formatter<'_>) -> fmt::Result {
let v = match self {
Self::None => "None",
Self::WaterOpc => "OPC water",
Self::WaterOpcSpecifyMolCount(c) => &format!("Water OPC. {c} mols"),
Self::Custom(_) => "Custom",
};
write!(f, "{v}")
}
}
impl PartialEq for Solvent {
fn eq(&self, other: &Self) -> bool {
match (self, other) {
(Self::WaterOpc, Self::WaterOpc) => true,
(Self::WaterOpcSpecifyMolCount(a), Self::WaterOpcSpecifyMolCount(b)) => a == b,
(Self::Custom((_, water_a)), Self::Custom((_, water_b))) => water_a == water_b,
_ => false,
}
}
}
#[cfg(feature = "encode")]
impl bincode::Encode for Solvent {
fn encode<E: bincode::enc::Encoder>(
&self,
encoder: &mut E,
) -> Result<(), bincode::error::EncodeError> {
match self {
Self::None | Self::WaterOpc => {
0u32.encode(encoder)?;
}
Self::WaterOpcSpecifyMolCount(count) => {
1u32.encode(encoder)?;
count.encode(encoder)?;
}
Self::Custom(_) => {
0u32.encode(encoder)?;
}
}
Ok(())
}
}
#[cfg(feature = "encode")]
impl<Context> bincode::Decode<Context> for Solvent {
fn decode<D: bincode::de::Decoder<Context = Context>>(
decoder: &mut D,
) -> Result<Self, bincode::error::DecodeError> {
let variant = u32::decode(decoder)?;
match variant {
0 => Ok(Self::WaterOpc),
1 => {
let count = usize::decode(decoder)?;
Ok(Self::WaterOpcSpecifyMolCount(count))
}
_ => Err(bincode::error::DecodeError::UnexpectedVariant {
type_name: "Solvent",
allowed: &bincode::error::AllowedEnumVariants::Range { min: 0, max: 1 },
found: variant,
}),
}
}
}
#[cfg(feature = "encode")]
impl<'de, Context> bincode::BorrowDecode<'de, Context> for Solvent {
fn borrow_decode<D: bincode::de::BorrowDecoder<'de, Context = Context>>(
decoder: &mut D,
) -> Result<Self, bincode::error::DecodeError> {
let variant = u32::borrow_decode(decoder)?;
match variant {
0 => Ok(Self::WaterOpc),
1 => {
let count = usize::borrow_decode(decoder)?;
Ok(Self::WaterOpcSpecifyMolCount(count))
}
_ => Err(bincode::error::DecodeError::UnexpectedVariant {
type_name: "Solvent",
allowed: &bincode::error::AllowedEnumVariants::Range { min: 0, max: 1 },
found: variant,
}),
}
}
}
#[derive(Copy, Clone, PartialEq)]
#[repr(u8)]
pub(crate) enum WaterSite {
O = 1,
M = 2,
H0 = 3,
H1 = 4,
}
#[derive(Clone, Copy, Default)]
pub struct ForcesOnWaterMol {
pub f_o: Vec3,
pub f_h0: Vec3,
pub f_h1: Vec3,
pub f_m: Vec3,
}
#[allow(unused)]
#[cfg(target_arch = "x86_64")]
#[derive(Clone, Copy, Default)]
pub struct ForcesOnWaterMolx8 {
pub f_o: Vec3x8,
pub f_h0: Vec3x8,
pub f_h1: Vec3x8,
pub f_m: Vec3x8,
}
#[allow(unused)]
#[cfg(target_arch = "x86_64")]
#[derive(Clone, Copy, Default)]
pub struct ForcesOnWaterMolx16 {
pub f_o: Vec3x16,
pub f_h0: Vec3x16,
pub f_h1: Vec3x16,
pub f_m: Vec3x16,
}
#[derive(Clone, Debug)]
pub struct WaterMol {
pub o: AtomDynamics,
pub h0: AtomDynamics,
pub h1: AtomDynamics,
pub m: AtomDynamics,
}
#[allow(unused)]
#[cfg(target_arch = "x86_64")]
pub struct WaterMolx8 {
pub o: AtomDynamicsx8,
pub h0: AtomDynamicsx8,
pub h1: AtomDynamicsx8,
pub m: AtomDynamicsx8,
}
#[allow(unused)]
#[cfg(target_arch = "x86_64")]
pub struct WaterMolx16 {
pub o: AtomDynamicsx16,
pub h0: AtomDynamicsx16,
pub h1: AtomDynamicsx16,
pub m: AtomDynamicsx16,
}
impl WaterMol {
pub fn new(o_pos: Vec3F32, vel: Vec3F32, orientation: QuaternionF32) -> Self {
let z_local = orientation.rotate_vec(Z_VEC);
let e_local = orientation.rotate_vec(X_VEC);
let h0_dir = (z_local * H_O_H_θ_HALF.cos() + e_local * H_O_H_θ_HALF.sin()).to_normalized();
let h1_dir = (z_local * H_O_H_θ_HALF.cos() - e_local * H_O_H_θ_HALF.sin()).to_normalized();
let h0_pos = o_pos + h0_dir * O_H_R;
let h1_pos = o_pos + h1_dir * O_H_R;
let ep_pos = o_pos + (h0_pos - o_pos + h1_pos - o_pos).to_normalized() * O_EP_R;
let h0 = AtomDynamics {
force_field_type: String::from("HW"),
element: Element::Hydrogen,
posit: h0_pos,
vel,
mass: H_MASS,
partial_charge: Q_H,
..Default::default()
};
Self {
o: AtomDynamics {
force_field_type: String::from("OW"),
posit: o_pos,
element: Element::Oxygen,
mass: O_MASS,
partial_charge: 0.,
lj_sigma: O_SIGMA,
lj_eps: O_EPS,
..h0.clone()
},
h1: AtomDynamics {
posit: h1_pos,
..h0.clone()
},
m: AtomDynamics {
force_field_type: String::from("EP"),
posit: ep_pos,
element: Element::Potassium, mass: 0.,
partial_charge: Q_EP,
..h0.clone()
},
h0,
}
}
pub(crate) fn _project_ep_force_to_real_sites(&mut self, cell: &SimBox) {
let r_O_H0 = self.o.posit + cell.min_image(self.h0.posit - self.o.posit) - self.o.posit;
let r_O_H1 = self.o.posit + cell.min_image(self.h1.posit - self.o.posit) - self.o.posit;
let s = r_O_H0 + r_O_H1;
let s_norm = s.magnitude();
if s_norm < 1e-12 {
self.o.force += self.m.force;
self.m.force = Vec3F32::new_zero();
return;
}
let f_m = self.m.force;
let u = s / s_norm;
let fm_parallel = u * f_m.dot(u);
let fm_perp = f_m - fm_parallel; let scale = O_EP_R / s_norm;
let fh = fm_perp * scale; let fo = f_m - fh * 2.0;
self.m.force = Vec3F32::new_zero();
self.o.force += fo;
self.h0.force += fh;
self.h1.force += fh;
}
pub(crate) fn project_ep_force_optimized(&mut self) {
let f_m = self.m.force;
self.o.force += f_m * C_O;
self.h0.force += f_m * C_H;
self.h1.force += f_m * C_H;
self.m.force = Vec3F32::new_zero();
}
pub(crate) fn update_virtual_site(&mut self) {
let v_h0 = self.h0.posit - self.o.posit;
let v_h1 = self.h1.posit - self.o.posit;
let bis = v_h0 + v_h1;
self.m.posit = self.o.posit + bis.to_normalized() * O_EP_R;
self.m.vel = (self.h0.vel + self.h1.vel) * 0.5;
}
}