YetAnotherSimulationSuite.CalculatorType
Calculator

Structure for a calculator object for MD.

Fields

  • b: Potential builder.
  • e: Energy function.
  • f!: Force function.
  • ef!: Energy and force function.
  • energy_unit: Unitful energy unit
  • force_unit: Unitful force unit
  • time_unit: Unitful time unit
  • constraints: Constraints.
source
YetAnotherSimulationSuite.CalculatorMethod
Calculator(b; E=nothing, F=nothing, EF=nothing, constraints=nothing)

Construct a Calculator object.

Arguments

  • b: Potential builder.
  • E: Energy function (optional).
  • F: Force function (optional).
  • EF: Energy and force function (optional).
  • constraints: Constraints (optional).

Returns

  • Calculator object.
source
YetAnotherSimulationSuite.CellType
Cell{D, B, I, F, S}

Structure representing a simulation cell.

Fields

  • lattice: Lattice matrix.
  • scaled_pos: Scaled positions.
  • velocity: Velocities.
  • masses: Masses.
  • symbols: Atomic symbols.
  • mask: Mask vector.
  • PBC: Periodic boundary conditions.
  • NC: Neighbor counts.
source
YetAnotherSimulationSuite.CellMethod
Cell(lat, spos, vel, mas, sym, mask, PBC, NC)

Construct a Cell object from lattice, positions, velocities, etc.

Arguments

  • lat: Lattice matrix.
  • spos: Scaled positions.
  • vel: Velocities.
  • mas: Masses.
  • sym: Symbols.
  • mask: Mask vector.
  • PBC: Periodic boundary conditions.
  • NC: Neighbor counts.

Returns

  • Cell object.
source
YetAnotherSimulationSuite.DynamicsType
Dynamics{T,D,B,P,PV,I,F,S}

Structure holding all MD simulation variables.

Fields

  • m: Masses.
  • s: Symbols.
  • mols: Molecule indices.
  • temp: Temperatures.
  • energy: Energies.
  • forces: Forces.
  • potVars: Potential variables.
  • PBC: Periodic boundary conditions.
  • NC: Neighbor counts.
  • ensemble: Ensemble object.
source
YetAnotherSimulationSuite.HiddenOptVarsType
HiddenOptVars

Structure for hidden variable optimization.

Fields

  • potVars: Potential variables.
  • cellBuf: Cell buffer.
  • superBuf: Supercell buffer.
  • scaleEnergy: Energy scaling factor.
  • PBC: Periodic boundary conditions.
  • mols: Molecule indices.
  • pars: Pair indices.
  • T: Transformation matrix.
source
YetAnotherSimulationSuite.ParticleMethod
Particle(r, v, m, s)

Construct a Particle from position, velocity, mass, and symbol.

Arguments

  • r: Position vector.
  • v: Velocity vector.
  • m: Mass.
  • s: Symbol.

Returns

  • Particle object.
source
YetAnotherSimulationSuite.TrajMethod
Traj(imgs, mas, sym, lat)

Construct a Traj object from images, masses, symbols, and lattice.

Arguments

  • imgs: Vector of Image objects.
  • mas: Vector of masses.
  • sym: Vector of symbols.
  • lat: Lattice matrix.

Returns

  • Traj object.
source
YetAnotherSimulationSuite.optVarsType
optVars

Structure holding optimization variables for geometry optimization.

Fields

  • potVars: Potential variables.
  • mols: Molecule indices.
  • m: Masses.
  • PBC: Periodic boundary conditions.
  • NC: Neighbor counts.
  • lattice: Lattice matrix.
source
YetAnotherSimulationSuite.vacfInpsType
vacfInps

Input structure for velocity autocorrelation function (VACF) calculations.

Fields

  • vel: Velocity data.
  • mas: Masses.
  • Hz: Sampling frequency.
  • norm: Normalize flag.
  • win: Window function.
  • pad: Padding factor.
  • mir: Mirror flag.
source
Base.deleteat!Method
deleteat!(cell::MyCell, iter)

Remove atoms at given indices from a cell.

Arguments

  • cell: MyCell object.
  • iter: Indices to remove.

Side Effects

  • Modifies the cell in-place.
source
Base.getindexMethod
getindex(cell::MyCell, inds)

Get cell atoms at inds indicies

Arguments

  • cell: MyCell object
  • inds: Indices to get

Returns

  • Cell object
source
Base.lengthMethod
length(cell::MyCell)

Get number of atoms in cell

Arguments

  • cell: MyCell object

Returns

  • Number of atoms in cell
source
Base.lengthMethod
length(traj::MyTraj)

Get the number of images in a trajectory.

Arguments

  • traj: Traj object.

Returns

  • Number of images.
source
Base.push!Method
push!(tj::MyTraj, solu::SciMLBase.ODESolution; dt=fs, step=1)

Append images from an ODE solution to a trajectory.

Arguments

  • tj: Traj object.
  • solu: ODE solution object.
  • dt: Time step (default: fs).
  • step: Step interval (default: 1).

Side Effects

  • Modifies tj in-place.
source
Base.repeatMethod
repeat(cell::MyCell, count::Integer)

Repeat a cell multiple times.

Arguments

  • cell: MyCell object.
  • count: Number of repetitions.

Returns

  • New repeated Cell object.
source
Base.runMethod
run(calc, cell, tmax, dt, ensemble; algo=VelocityVerlet(), split=1, kwargs...)

Run an MD simulation for a cell.

Arguments

  • calc: Calculator object.
  • cell: MyCell object.
  • tmax: Time span tuple.
  • dt: Time step.
  • ensemble: Ensemble object.
  • algo: ODE solver algorithm (default: VelocityVerlet()).
  • split: Number of segments (default: 1).
  • kwargs: Additional keyword arguments.

Returns

  • Processed trajectory or solution.
source
Base.runMethod
run(calc, cell, tspan, dt, ensemble; algo=VelocityVerlet(), split=1, kwargs...)

Run an MD simulation for a cell.

Arguments

  • calc: Calculator object.
  • cell: MyCell object.
  • tspan: Time span tuple.
  • dt: Time step.
  • ensemble: Ensemble object.
  • algo: ODE solver algorithm (default: VelocityVerlet()).
  • split: Number of segments (default: 1).
  • kwargs: Additional keyword arguments.

Returns

  • Processed trajectory or solution.
source
Base.runMethod
run(calc, bdys, tmax, dt, ensemble; algo=VelocityVerlet(), split=1, kwargs...)

Run an MD simulation for a set of atoms.

Arguments

  • calc: Calculator object.
  • bdys: Vector of MyAtoms.
  • tmax: Time length of simulation.
  • dt: Time step.
  • ensemble: Ensemble object.
  • algo: ODE solver algorithm (default: VelocityVerlet()).
  • split: Number of segments (default: 1).
  • kwargs: Additional keyword arguments.

Returns

  • Processed trajectory or solution.
source
Base.runMethod
run(calc, bdys, tspan, dt, ensemble; algo=VelocityVerlet(), split=1, kwargs...)

Run an MD simulation for a set of atoms.

Arguments

  • calc: Calculator object.
  • bdys: Vector of MyAtoms.
  • tspan: Time span tuple.
  • dt: Time step.
  • ensemble: Ensemble object.
  • algo: ODE solver algorithm (default: VelocityVerlet()).
  • split: Number of segments (default: 1).
  • kwargs: Additional keyword arguments.

Returns

  • Processed trajectory or solution.
source
Base.showMethod
show(io::IO, cell::MyCell)

Display cell struct information

Arguments

  • io: IO
  • cell: MyCell object
source
Base.showMethod
show(io::IO, traj::MyTraj)

Custom display for MyTraj objects.

Arguments

  • io: IO stream.
  • traj: MyTraj object.
source
Base.writeMethod
Base.write(file::String, cell::MyCell)

Write a MyCell object to a file using Chemfiles.

Arguments

  • file: Output file path.
  • cell: MyCell object to write.

Side Effects

  • Writes atomic coordinates, velocities, and cell information to file.
source
Base.writeMethod
Base.write(file::String, traj::MyTraj; step=1)

Write a trajectory to a file using Chemfiles.

Arguments

  • file: Output file path.
  • traj: MyTraj object to write.
  • step: Step interval for writing frames (default: 1).

Side Effects

  • Writes trajectory frames, including positions, velocities, energies, and forces, to file.
source
Base.writeMethod
Base.write(file::String, bdys::Vector{MyAtoms})

Write a vector of MyAtoms to a file using Chemfiles.

Arguments

  • file: Output file path.
  • bdys: Vector of MyAtoms to write.

Side Effects

  • Writes atomic coordinates and velocities to file.
source
LinearAlgebra.rotate!Method
LinearAlgebra.rotate!(bdys::Vector{MyAtoms}, R::Matrix, about::Vector{Float64})

Rotate all atoms in bdys by rotation matrix R about a point about in-place.

Arguments

  • bdys: Vector of MyAtoms objects.
  • R: 3x3 rotation matrix.
  • about: Point to rotate about.
source
LinearAlgebra.rotate!Method
LinearAlgebra.rotate!(bdys::Vector{MyAtoms}, R::Matrix)

Rotate all atoms in bdys by rotation matrix R in-place.

Arguments

  • bdys: Vector of MyAtoms objects.
  • R: 3x3 rotation matrix.
source
LinearAlgebra.rotate!Method
LinearAlgebra.rotate!(bdys::Vector{MyAtoms}, abc::Tuple{F,F,F}, about::Vector{F}) where F<:AbstractFloat

Rotate all atoms in bdys by Euler angles (α, β, γ) about a point about in-place.

Arguments

  • bdys: Vector of MyAtoms objects.
  • abc: Tuple of Euler angles (α, β, γ) in radians.
  • about: Point to rotate about.
source
LinearAlgebra.rotate!Method
LinearAlgebra.rotate!(bdys::Vector{MyAtoms}, abc::Tuple{F,F,F}) where F <: AbstractFloat

Rotate all atoms in bdys by Euler angles (α, β, γ) in-place.

Arguments

  • bdys: Vector of MyAtoms objects.
  • abc: Tuple of Euler angles (α, β, γ) in radians.
source
YetAnotherSimulationSuite.CoMMethod
CoM(pos, mas)

Compute the center of mass from positions and masses.

Arguments

  • pos: Vector of position vectors.
  • mas: Vector of masses.

Returns

  • Center of mass as a vector.
source
YetAnotherSimulationSuite.CoMMethod
CoM(bdys)

Compute the center of mass for a collection of atoms.

Arguments

  • bdys: Vector of objects with fields m (mass) and r (position).

Returns

  • Center of mass as a 3-element vector.
source
YetAnotherSimulationSuite.HGNNMethod

CO-CO Potential

Based on: Chen, Jun, et al. "Energy transfer between vibrationally excited carbon monoxide based on a highly accurate six-dimensional potential energy surface." The Journal of Chemical Physics 153.5 (2020).

Link: https://pubs.aip.org/aip/jcp/article/153/5/054310/1065758

source
YetAnotherSimulationSuite.MBXMethod

MBX potentials

This allows uing all the potentials found within the MBX package.

MBX GitHub: https://github.com/paesanilab/MBX/tree/master

MBX Paper: https://pubs.aip.org/aip/jcp/article/159/5/054802/2904909/MBX-A-many-body-energy-and-force-calculator-for

source
YetAnotherSimulationSuite.MvHffMethod

CO-CO Potential

Based on: van Hemert, Marc C., Junko Takahashi, and Ewine F. van Dishoeck. "Molecular dynamics study of the photodesorption of CO ice." The Journal of Physical Chemistry A 119.24 (2015): 6354-6369.

Link: https://pubs.acs.org/doi/full/10.1021/acs.jpca.5b02611

source
YetAnotherSimulationSuite.SCMEfMethod

SCME/f potential

This calls an ASE calculator of the SCME/f potential (not great) I still want to make a better wrapper. Idealy, one that calls the c++ routines directly (like MBX).

SCME Paper: https://pubs.rsc.org/en/content/articlehtml/2013/cp/c3cp52097h

SCME/f Paper: https://pubs.acs.org/doi/full/10.1021/acs.jctc.2c00598

source
YetAnotherSimulationSuite.SPCMethod

Simple Point Charge (SPC) Models

SPC/F

Based on: Toukan, Kahled, and Aneesur Rahman. "Molecular-dynamics study of atomic motions in water." Physical Review B 31.5 (1985): 2643.

Link: https://journals.aps.org/prb/abstract/10.1103/PhysRevB.31.2643

SPC/Fd

Based on: Dang, Liem X., and B. Montgomery Pettitt. "Simple intramolecular model potentials for water." Journal of physical chemistry 91.12 (1987): 3349-3354.

Link: https://doi.org/10.1021/j100296a048

SPC/Fw

Based on: Wu, Yujie, Harald L. Tepper, and Gregory A. Voth. "Flexible simple point-charge water model with improved liquid-state properties." The Journal of chemical physics 124.2 (2006).

Link: https://doi.org/10.1063/1.2136877

source
YetAnotherSimulationSuite.TIP4PfMethod

TIP4P/2005f

Based on: González, M. A., & Abascal, J. L. (2011). A flexible model for water based on TIP4P/2005. The Journal of chemical physics, 135(22).

Link: https://pubs.aip.org/aip/jcp/article/135/22/224516/190786

source
YetAnotherSimulationSuite.VDOSMethod
VDOS(inp; atms=nothing)

Compute the vibrational density of states (VDOS) from VACF input.

Arguments

  • inp: vacfInps object.
  • atms: (Optional) Indices of atoms to include.

Returns

  • vacfOut object with VDOS data.
source
YetAnotherSimulationSuite.adfMethod

Compute the angular distribution function (ADF) between molecular orientations.

Arguments

  • bdys::Vector{MyAtoms}: Vector of MyAtoms objects.
  • kwargs...: Additional keyword arguments for kde_lscv.

Returns

  • (x, y): Tuple of angle values (in degrees) and ADF.
source
YetAnotherSimulationSuite.animateModeMethod
animateMode(bdys::Vector{MyAtoms}, mode, fileName; c=1.0)

Animate a vibrational mode and write the trajectory to a file.

Arguments

  • bdys: Vector of MyAtoms objects.
  • mode: Mode vector to animate.
  • fileName: Output file name.
  • c: (Optional) Amplitude scaling factor (default: 1.0).
source
YetAnotherSimulationSuite.centerMethod
center(bdys::Vector{MyAtoms})

Compute the geometric center of a set of atoms.

Arguments

  • bdys: Vector of MyAtoms.

Returns

  • Center position as a vector.
source
YetAnotherSimulationSuite.countNearestNeighborsMethod
countNearestNeighbors(mol, bdys; rmax=4.5)

Count the number of neighboring molecules within a cutoff distance.

Arguments

  • mol: The molecule of interest.
  • bdys: Collection of molecules to search for neighbors.
  • rmax: (Optional) Cutoff distance for neighbor search (default 4.5).

Returns

  • Number of neighboring molecules within rmax of mol.
source
YetAnotherSimulationSuite.densityMethod

Compute the average number density of atoms in a given set.

Arguments

  • bdys::Vector{MyAtoms}: Vector of MyAtoms objects representing atoms.

Returns

  • Float64: Average number density (atoms per unit volume) within the maximum radial extent from the center.
source
YetAnotherSimulationSuite.densityMethod

Compute the average number density for a set of positions relative to a given center.

Arguments

  • o::Vector{Float64}: The center point as a vector.
  • pts::Vector{A}: Vector of position vectors (any subtype of AbstractArray).

Returns

  • Float64: Average number density (points per unit volume) within the maximum radial extent from the center.
source
YetAnotherSimulationSuite.diffDotSqrtMethod
diffDotSqrt(v2, v1)

Compute the vector difference and its Euclidean norm between two vectors.

Arguments

  • v2, v1: Input vectors.

Returns

  • Tuple: (Euclidean distance, difference vector).
source
YetAnotherSimulationSuite.doRunMethod
doRun(calc, vel, pos, tspan, simu, algo, dt, split; kwargs...)

Run an MD simulation, optionally splitting into segments.

Arguments

  • calc: Calculator object.
  • vel: Initial velocities.
  • pos: Initial positions.
  • tspan: Time span tuple.
  • simu: Dynamics object.
  • algo: ODE solver algorithm.
  • dt: Time step.
  • split: Number of segments.
  • kwargs: Additional keyword arguments.

Returns

  • Processed trajectory or solution.
source
YetAnotherSimulationSuite.dyn!Method
dyn!(dv, v, u, p, t, calc)

Compute the time derivative for MD integration.

Arguments

  • dv: Output derivative.
  • v: Velocities.
  • u: Positions.
  • p: Dynamics object.
  • t: Time.
  • calc: Calculator object.

Side Effects

  • Modifies dv in-place.
source
YetAnotherSimulationSuite.fg!Method
fg!(F, G, x, p, calc::MyCalc)

Evaluate energy and/or forces and gradients for optimization.

Arguments

  • F: If not nothing, return energy.
  • G: If not nothing, fill with gradients.
  • x: Coordinate vector.
  • p: Potential variables.
  • calc: Calculator object.

Returns

  • Energy if F is not nothing.
source
YetAnotherSimulationSuite.findPeaksMethod
findPeaks(arr; min=0.0, max=1e30, width=1)

Find indices of local maxima (peaks) in an array.

Arguments

  • arr: Array of numeric values to search for peaks.
  • min: Minimum value threshold for a peak (default: 0.0).
  • max: Maximum value threshold for a peak (default: 1e30).
  • width: Minimum width for a peak (default: 1).

Returns

  • Vector of indices where peaks are found.
source
YetAnotherSimulationSuite.findTurningPointsMethod
findTurningPoints(arr)

Find indices of local maxima and minima (turning points) in an array.

Arguments

  • arr: Array of numeric values.

Returns

  • Vector of indices where turning points are found.
source
YetAnotherSimulationSuite.getDistanceMatrixMethod
getDistanceMatrix(cell::MyCell)

Compute the pairwise distance matrix for atoms in a cell, accounting for periodicity.

Arguments

  • cell: MyCell object.

Returns

  • Symmetric matrix of distances between atoms.
source
YetAnotherSimulationSuite.getForcesMethod
getForces(calc::MyCalc, cell::MyCell)

Compute the forces on all atoms in a cell using a calculator.

Arguments

  • calc: Calculator object.
  • cell: MyCell object.

Returns

  • Vector of force vectors.
source
YetAnotherSimulationSuite.getForcesMethod
getForces(calc::MyCalc, bdys::Vector{MyAtoms})

Compute the forces on a set of atoms using a calculator.

Arguments

  • calc: Calculator object.
  • bdys: Vector of MyAtoms objects.

Returns

  • Vector of force vectors.
source
YetAnotherSimulationSuite.getHarmonicFreqsMethod
getHarmonicFreqs(calc::MyCalc, obj::Union{MyCell, Vector{MyAtoms}})

Compute harmonic vibrational frequencies and modes for a cell or molecule.

Arguments

  • calc: Calculator object (MyCalc).
  • obj: MyCell or vector of MyAtoms.

Returns

  • Tuple: (vector of frequencies, matrix of modes).
source
YetAnotherSimulationSuite.getIPRMethod
getIPR(modes; N=2)

Compute the inverse participation ratio (IPR) for vibrational modes.

Arguments

  • modes: Matrix of mode vectors.
  • N: Number of atoms per molecule (default: 2).

Returns

  • Matrix of IPR values.
source
YetAnotherSimulationSuite.getImageMethod
getImage(solu::SciMLBase.ODESolution, i::Int, dt::Float64)

Extract an Image from an ODE solution at a given index.

Arguments

  • solu: ODE solution object.
  • i: Index of the time step.
  • dt: Time step size.

Returns

  • Image object.
source
YetAnotherSimulationSuite.getModeInteractionPESMethod
getModeInteractionPES(EoM!, bdys, mode; range=collect(-1:0.001:2))

Compute the interaction potential energy surface along a vibrational mode.

Arguments

  • EoM!: Energy function.
  • bdys: Vector of MyAtoms objects.
  • mode: Mode vector.
  • range: (Optional) Range of displacements (default: -1:0.001:2).

Returns

  • Tuple: (displacement values, interaction energy values).
source
YetAnotherSimulationSuite.getModePESMethod
getModePES(EoM, bdys, mode; range=collect(-1:0.001:2))

Compute the potential energy surface (PES) along a vibrational mode.

Arguments

  • EoM: Energy function.
  • bdys: Vector of MyAtoms objects.
  • mode: Mode vector.
  • range: (Optional) Range of displacements (default: -1:0.001:2).

Returns

  • Tuple: (displacement values, energy values).
source
YetAnotherSimulationSuite.getMolsMethod
getMols(cell::MyCell, rmax::Float64)

Cluster atoms in a cell into molecules using DBSCAN with cutoff rmax.

Arguments

  • cell: MyCell object.
  • rmax: Distance cutoff for clustering.

Returns

  • Vector of vectors, each containing indices of atoms in a molecule.
source
YetAnotherSimulationSuite.getMolsMethod
getMols(bdys::Vector{MyAtoms}, rmax::Float64)

Cluster atoms in bdys into molecules using DBSCAN with cutoff rmax.

Arguments

  • bdys: Vector of MyAtoms objects.
  • rmax: Distance cutoff for clustering.

Returns

  • Vector of vectors, each containing indices of atoms in a molecule.
source
YetAnotherSimulationSuite.getNewBdysMethod
getNewBdys(bdys, res)

Construct new atom objects from optimization results.

Arguments

  • bdys: Original vector of MyAtoms.
  • res: Optimization result.

Returns

  • Vector of new MyAtoms objects.
source
YetAnotherSimulationSuite.getNumericalStressMethod
getNumericalStress(calc::MyCalc, cell::MyCell; eps=1e-6)

Compute the numerical stress tensor for a cell using finite differences.

Arguments

  • calc: Calculator object (MyCalc).
  • cell: Cell object (MyCell).
  • eps: Strain increment for finite difference (default: 1e-6).

Returns

  • 3x3 stress tensor matrix.
source
YetAnotherSimulationSuite.getNumericalStressOrthogonalMethod
getNumericalStressOrthogonal(calc::MyCalc, cell::MyCell; eps=1e-6)

Compute the numerical stress tensor for a cell using finite differences, only for diagonal terms (orthogonal cell).

Arguments

  • calc: Calculator object (MyCalc).
  • cell: Cell object (MyCell).
  • eps: Strain increment for finite difference (default: 1e-6).

Returns

  • 3x3 stress tensor matrix (only diagonal terms are nonzero).
source
YetAnotherSimulationSuite.getPRMethod
getPR(ipr; x=1.6e-5)

Compute the participation ratio (PR) from IPR data.

Arguments

  • ipr: IPR matrix.
  • x: Threshold value (default: 1.6e-5).

Returns

  • Vector of participation ratios.
source
YetAnotherSimulationSuite.getPotEnergyMethod
getPotEnergy(calc::MyCalc, obj::Union{MyCell, Vector{MyAtoms}})

Compute the potential energy of a cell or molecule.

Arguments

  • calc: Calculator object (MyCalc).
  • obj: MyCell or vector of MyAtoms.

Returns

  • Potential energy (Float64).
source
YetAnotherSimulationSuite.getRdfMethod

Compute the radial distribution function (RDF) from a list of distances and a density.

Arguments

  • d: Array of interatomic distances.
  • ρ: Number density.
  • kwargs...: Additional keyword arguments passed to kde_lscv.

Returns

  • (x, y): Tuple of radius values x and normalized RDF values y.
source
YetAnotherSimulationSuite.getRotationMatrixMethod
getRotationMatrix(α::F, β::F, γ::F) where F <: AbstractFloat

Construct a rotation matrix from Euler angles (α, β, γ).

Arguments

  • α, β, γ: Euler angles in radians.

Returns

  • 3x3 rotation matrix.
source
YetAnotherSimulationSuite.getScaledPos!Method
getScaledPos!(cell, x0)

Update scaled positions in a cell from Cartesian coordinates.

Arguments

  • cell: MyCell object.
  • x0: Cartesian coordinates.

Side Effects

  • Modifies the cell in-place.
source
YetAnotherSimulationSuite.getScaledPosMethod
getScaledPos(x0, lattice)

Get scaled positions from Cartesian coordinates and lattice.

Arguments

  • x0: Cartesian coordinates.
  • lattice: Lattice matrix.

Returns

  • Vector of scaled positions.
source
YetAnotherSimulationSuite.getScaledPosMethod
getScaledPos(bdys::Vector{MyAtoms}, lattice)

Get scaled positions for a set of atoms given a lattice.

Arguments

  • bdys: Vector of MyAtoms.
  • lattice: Lattice matrix.

Returns

  • Vector of scaled positions.
source
YetAnotherSimulationSuite.getSurfaceMoleculesMethod
getSurfaceMolecules(bdys::Vector{MyAtoms}; α=nothing)

Get the surface molecules from a set of atoms using alpha shapes.

Arguments

  • bdys: Vector of MyAtoms.
  • α: Alpha parameter (optional).

Returns

  • Vector of surface MyAtoms.
source
YetAnotherSimulationSuite.getVibEnergyMethod
getVibEnergy(mol::Vector{MyAtoms}, eignvec; calc=nothing)

Compute the vibrational energy of a molecule along a mode.

Arguments

  • mol: Vector of MyAtoms.
  • eignvec: Mode vector.
  • calc: (Optional) Calculator for potential energy.

Returns

  • Vibrational energy (Float64).
source
YetAnotherSimulationSuite.hiddenEoMMethod
hiddenEoM(F, G, Γ, calc::MyCalc, vars, x)

Energy and gradient evaluation for hidden variable optimization.

Arguments

  • F: Energy output (or nothing).
  • G: Gradient output (or nothing).
  • Γ: Gradient buffer.
  • calc: Calculator object.
  • vars: HiddenOptVars object.
  • x: Coordinate vector.

Returns

  • Energy value if F is not nothing.
source
YetAnotherSimulationSuite.hiddenOptMethod
hiddenOpt(calc::MyCalc, algo, cell, T; kwargs...)

Optimize a cell using a hidden variable approach and a supercell transformation.

Arguments

  • calc: Calculator object.
  • algo: Optimization algorithm.
  • cell: MyCell object.
  • T: Transformation matrix.
  • kwargs: Additional options.

Returns

  • Optimized MyCell object.
source
YetAnotherSimulationSuite.ifPropertyMethod
ifProperty(frame, props, p)

Return the value of property p from frame if it exists in props, else return 0.0.

Arguments

  • frame: Chemfiles Frame object.
  • props: List of property names.
  • p: Property name to look for.

Returns

  • Value of the property as Float64, or 0.0 if not found.
source
YetAnotherSimulationSuite.makeBdysMethod
makeBdys(tj::MyTraj, i::Int)

Construct a vector of MyAtoms from a trajectory at a given image index.

Arguments

  • tj: Traj object.
  • i: Image index.

Returns

  • Vector of MyAtoms.
source
YetAnotherSimulationSuite.makeCellMethod
makeCell(bdys, lattice; mask, PBC, NC)

Construct a Cell from a vector of MyAtoms and a lattice.

Arguments

  • bdys: Vector of MyAtoms.
  • lattice: Lattice matrix.
  • mask: Mask vector (optional).
  • PBC: Periodic boundary conditions (optional).
  • NC: Neighbor counts (optional).

Returns

  • Cell object.
source
YetAnotherSimulationSuite.makeSuperCell!Method
makeSuperCell!(super, cell, T)

In-place creation of a supercell from a cell and transformation matrix.

Arguments

  • super: Supercell object to fill.
  • cell: Cell to replicate.
  • T: Transformation matrix.

Side Effects

  • Modifies the supercell in-place.
source
YetAnotherSimulationSuite.mkvarMethod
mkvar(x)

Evaluate a symbol from a string and return its value.

Arguments

  • x: String representing the variable name.

Returns

  • Value of the variable with name x.
source
YetAnotherSimulationSuite.myRepeatMethod
myRepeat(A::Vector{MVector{D,Float64}}, count::Integer, mask::Vector{Bool}) where D

Repeat elements of a vector of static vectors, excluding masked elements, and return a deep copy.

Arguments

  • A: Vector of MVector{D,Float64}.
  • count: Number of times to repeat unmasked elements.
  • mask: Boolean mask vector (true = keep, false = repeat).

Returns

  • Vector of repeated and copied static vectors.
source
YetAnotherSimulationSuite.myRepeatMethod
myRepeat(A::Vector{T}, count::Integer, mask::Vector{Bool}) where T <: Union{String, Number}

Repeat elements of a vector, excluding masked elements, and concatenate with the original.

Arguments

  • A: Input vector.
  • count: Number of times to repeat unmasked elements.
  • mask: Boolean mask vector (true = keep, false = repeat).

Returns

  • Concatenated vector with repeated elements.
source
YetAnotherSimulationSuite.optMethod
opt(calc::MyCalc, algo, cell::MyCell; kwargs...)

Optimize the geometry of a cell.

Arguments

  • calc: Calculator object.
  • algo: Optimization algorithm.
  • cell: MyCell object.
  • kwargs: Additional options.

Returns

  • Optimized MyCell object.
source
YetAnotherSimulationSuite.optMethod
opt(calc::MyCalc, algo, bdys::Vector{MyAtoms}; kwargs...)

Optimize the geometry of a set of atoms.

Arguments

  • calc: Calculator object.
  • algo: Optimization algorithm.
  • bdys: Vector of MyAtoms objects.
  • kwargs: Additional options.

Returns

  • Optimized vector of MyAtoms.
source
YetAnotherSimulationSuite.optCellMethod
optCell(calc::MyCalc, algo, cell::MyCell; precon=nothing, kwargs...)

Optimize the lattice parameters of a cell.

Arguments

  • calc: Calculator object.
  • algo: Optimization algorithm.
  • cell: MyCell object.
  • precon: (Optional) Preconditioner.
  • kwargs: Additional options.

Returns

  • Optimized MyCell object.
source
YetAnotherSimulationSuite.periodicDistanceMethod
periodicDistance(x::AbstractArray, y::AbstractArray, vects)

Compute the minimum periodic distance between two points given lattice vectors.

Arguments

  • x, y: Position vectors.
  • vects: Iterable of lattice vectors.

Returns

  • Minimum periodic distance (Float64).
source
YetAnotherSimulationSuite.pickRandomMolMethod
pickRandomMol(bdys::Vector{MyAtoms}, loc)

Pick a random molecule from the surface or bulk.

Arguments

  • bdys: Vector of MyAtoms.
  • loc: "surf" or "bulk".

Returns

  • Vector of MyAtoms for the selected molecule.
source
YetAnotherSimulationSuite.prep4potMethod
prep4pot(builder, cell::MyCell)

Prepare variables for potential energy calculation from a cell.

Arguments

  • builder: Potential builder function.
  • cell: MyCell object.

Returns

  • Tuple: (coordinate vector, optVars object).
source
YetAnotherSimulationSuite.prep4potMethod
prep4pot(builder, bdys::Vector{MyAtoms})

Prepare variables for potential energy calculation from atoms.

Arguments

  • builder: Potential builder function.
  • bdys: Vector of MyAtoms objects.

Returns

  • Tuple: (coordinate vector, optVars object).
source
YetAnotherSimulationSuite.prepX0Method
prepX0(bdys::Vector{MyAtoms})

Prepare the initial coordinate vector from a set of atoms.

Arguments

  • bdys: Vector of MyAtoms objects.

Returns

  • Flat vector of coordinates.
source
YetAnotherSimulationSuite.processDynamicsMethod
processDynamics(solu::SciMLBase.ODESolution; dt=fs, step=1)

Convert an ODE solution to a Traj object.

Arguments

  • solu: ODE solution object.
  • dt: Time step (default: fs).
  • step: Step interval (default: 1).

Returns

  • Traj object.
source
YetAnotherSimulationSuite.processTmpFilesMethod
processTmpFiles(files; kwargs...)

Process a list of temporary files containing ODE solutions into a single trajectory.

Arguments

  • files: List of file paths.
  • kwargs: Additional keyword arguments.

Returns

  • Traj object.
source
YetAnotherSimulationSuite.randRotate!Method
randRotate!(bdys::Vector{MyAtoms}, about::Vector{Float64})

Apply a random rotation about a point about to all atoms in bdys in-place.

Arguments

  • bdys: Vector of MyAtoms objects.
  • about: Point to rotate about.
source
YetAnotherSimulationSuite.rdfMethod

Compute the radial distribution function for pairs of two specified atomic species.

Arguments

  • bdys::Vector{MyAtoms}: Vector of MyAtoms objects.
  • P::Pair{String, String}: Pair of atomic species labels (A, B).
  • kwargs...: Additional keyword arguments for getRdf.

Returns

  • (x, y): Tuple of radius values and RDF for the species pair.
source
YetAnotherSimulationSuite.rdfMethod

Compute the radial distribution function for all pairs of a given atomic species.

Arguments

  • bdys::Vector{MyAtoms}: Vector of MyAtoms objects.
  • A::String: Atomic species label.
  • kwargs...: Additional keyword arguments for getRdf.

Returns

  • (x, y): Tuple of radius values and RDF for the given species.
source
YetAnotherSimulationSuite.rdfMethod

Compute the radial distribution function between molecular centers of mass.

Arguments

  • bdys::Vector{MyAtoms}: Vector of MyAtoms objects.
  • rmol: Distance threshold for molecular grouping (default 1.2).
  • kwargs...: Additional keyword arguments for getRdf.

Returns

  • (x, y): Tuple of radius values and RDF.
source
YetAnotherSimulationSuite.readFrameMethod
readFrame(frame::Frame)

Convert a Chemfiles Frame to a vector of MyAtoms or a MyCell if lattice is present.

Arguments

  • frame: Chemfiles Frame object.

Returns

  • Vector of MyAtoms or MyCell object.
source
YetAnotherSimulationSuite.readImageMethod
readImage(frame::Frame)

Read an image (snapshot) from a Chemfiles Frame, extracting positions, velocities, forces, and properties.

Arguments

  • frame: Chemfiles Frame object.

Returns

  • Image object containing positions, velocities, time, temperature, energy, and forces.
source
YetAnotherSimulationSuite.readSystemMethod
readSystem(file::String)

Read a system from a file and return either a trajectory or a single frame.

Arguments

  • file: Path to the file.

Returns

  • Traj object if multiple frames, otherwise a single frame as MyCell or MyAtoms.
source
YetAnotherSimulationSuite.readTrajMethod
readTraj(buf::Trajectory, N::Int)

Read a trajectory from a Chemfiles Trajectory buffer.

Arguments

  • buf: Chemfiles Trajectory object.
  • N: Number of frames to read.

Returns

  • Traj object containing all images, masses, symbols, and lattice.
source
YetAnotherSimulationSuite.reorder!Method
reorder!(cell::MyCell, order::Vector{Int})

Reorder the atoms in a cell according to a given order.

Arguments

  • cell: MyCell object.
  • order: New order of indices.

Side Effects

  • Modifies the cell in-place.
source
YetAnotherSimulationSuite.replicate!Method
replicate!(super, cell, N)

Replicate a cell into a supercell.

Arguments

  • super: Supercell object to fill.
  • cell: Cell to replicate.
  • N: Replication factors along each axis.

Side Effects

  • Modifies the supercell in-place.
source
YetAnotherSimulationSuite.savGolMethod
savGol(y, ws, order; deriv=0)

Apply a Savitzky-Golay filter to the input signal y.

Arguments

  • y: Array of data points to be smoothed or differentiated.
  • ws: Window size (must be odd and ≥ 1).
  • order: Polynomial order for the filter (must satisfy ws ≥ order + 2).
  • deriv: (Optional) Order of the derivative to compute (default is 0, i.e., smoothing).

Returns

  • Filtered (or differentiated) signal as an array of the same length as y.

Throws

  • ArgumentError if input arguments are invalid.

Example

smoothed = savGol(data, 5, 2)
source
YetAnotherSimulationSuite.singleRunMethod
singleRun(calc, vel, pos, tspan, simu, algo, dt; kwargs...)

Run a single MD simulation segment.

Arguments

  • calc: Calculator object.
  • vel: Initial velocities.
  • pos: Initial positions.
  • tspan: Time span tuple.
  • simu: Dynamics object.
  • algo: ODE solver algorithm.
  • dt: Time step.
  • kwargs: Additional keyword arguments.

Returns

  • ODE solution object.
source
YetAnotherSimulationSuite.st!Method
st!(F, G, x, cell::MyCell, calc::MyCalc; precon=nothing, diag=false)

Evaluate stress and/or energy for cell optimization.

Arguments

  • F: If not nothing, return energy.
  • G: If not nothing, fill with gradients.
  • x: Lattice vector.
  • cell: MyCell object.
  • calc: Calculator object.
  • precon: Preconditioner (optional).
  • diag: Use only diagonal terms (default: false).

Returns

  • Energy if F is not nothing.
source
YetAnotherSimulationSuite.swapAtoms!Method
swapAtoms!(bdys::Vector{MyAtoms}, i, j)

Swap the positions of two atoms in a vector.

Arguments

  • bdys: Vector of MyAtoms.
  • i: Index of first atom.
  • j: Index of second atom.

Side Effects

  • Modifies the positions in-place.
source
YetAnotherSimulationSuite.swapIso!Method
swapIso!(bdys::Vector{MyAtoms}, swap, mas)

Swap isotopes (masses) for a set of atoms.

Arguments

  • bdys: Vector of MyAtoms.
  • swap: Indices to swap.
  • mas: New masses.

Side Effects

  • Modifies masses in-place.
source
YetAnotherSimulationSuite.transExcite!Method
transExcite!(mol::Vector{MyAtoms}, ke)

Excite the translational motion of a molecule to a given kinetic energy.

Arguments

  • mol: Vector of MyAtoms.
  • ke: Target kinetic energy.

Side Effects

  • Modifies velocities of atoms in-place.
source
YetAnotherSimulationSuite.translate!Method
translate!(bdys::Vector{MyAtoms}, v::Vector{Float64})

Translate all atoms in bdys by vector v in-place.

Arguments

  • bdys: Vector of MyAtoms objects.
  • v: Translation vector.
source
YetAnotherSimulationSuite.vCoMMethod
vCoM(vel, mas)

Compute the center of mass velocity from velocities and masses.

Arguments

  • vel: Vector of velocity vectors.
  • mas: Vector of masses.

Returns

  • Center of mass velocity as a vector.
source
YetAnotherSimulationSuite.vCoMMethod
vCoM(bdys)

Compute the center of mass velocity for a collection of atoms.

Arguments

  • bdys: Vector of objects with fields m (mass) and v (velocity).

Returns

  • Center of mass velocity as a vector.
source
YetAnotherSimulationSuite.vacf!Method
vacf!(inp, out; atms=nothing)

Compute the velocity autocorrelation function (VACF).

Arguments

  • inp: vacfInps object.
  • out: vacfOut object.
  • atms: (Optional) Indices of atoms to include.
source
YetAnotherSimulationSuite.vibExcite!Method
vibExcite!(mol::Vector{MyAtoms}, eignvec, E)

Excite a vibrational mode of a molecule to a given energy.

Arguments

  • mol: Vector of MyAtoms.
  • eignvec: Mode vector.
  • E: Target energy.

Side Effects

  • Modifies velocities of atoms in-place.
source
YetAnotherSimulationSuite.wrap!Method
wrap!(bdys::Vector{MyAtoms}, lattice)

Wrap all atoms into the primary unit cell.

Arguments

  • bdys: Vector of MyAtoms.
  • lattice: Lattice matrix.

Side Effects

  • Modifies positions in-place.
source
YetAnotherSimulationSuite.zeroVCoM!Method
zeroVCoM!(bdys)

Remove the center of mass velocity from a collection of atoms in-place.

Arguments

  • bdys: Vector of objects with fields m (mass) and v (velocity).

Side Effects

  • Modifies velocities in-place to set total momentum to zero.
source