API Reference
This page provides a list of all documented types and functions and in Atomistic.jl.
If you are looking for more specifics on the Atomistic API, see the Implementing the Atomistic API page.
Atomistic.MolecularDynamicsResult — TypeMolecularDynamicsResultAbstract type to be extended by all concrete structs representing the result of a molecular dynamics simulation.
Atomistic.MolecularDynamicsSimulator — TypeMolecularDynamicsSimulatorAbstract type to be extended by all concrete molecular dynamics simulators.
Atomistic.MollyResult — TypeMollyResult <: MolecularDynamicsResultThe result generated from running an MollySimulator.
Field descriptions
system::SystemtheSystemresulting from the simulationsimulator::MollySimulatortheMollySimulatorthat produced the result
Atomistic.MollySimulator — TypeMollySimulator{S,T<:Unitful.Time,C} <: MolecularDynamicsSimulatorA wrapper around Molly to implement the Atomistic API.
Type parameter descriptions
Sthe type of simulatorT<:Unitful.Timethe type for the time fieldsCthe type for the coupling field
Field descriptions
Δt::Tthe time between timesteps, assumed to be in atomic unitssteps::Integerthe number of timesteps for the simulationt₀::Tthe starting time of the simulation, assumed to be in atomic units; defaults to 0coupling::Cthe coupling for the simulation; many options are defined byMolly, but a user could also define a custom thermostat; defaults toNoCoupling()stride::Intthe number of timesteps between logger runs defaults to 1
Atomistic.NBSResult — TypeNBSResult <: MolecularDynamicsResultThe result generated from running an NBSimulator.
Field descriptions
result::SimulationResultthe standard simulation result fromNBodySimulatorenergy_cache::Vector{Float64}cached potential energy values from the simulation
Atomistic.NBSimulator — TypeNBSimulator <: MolecularDynamicsSimulatorA wrapper around NBodySimulator to implement the Atomistic API.
Type parameter descriptions
Sthe type of ODE solverT<:Unitful.Timethe type for the time fieldsCthe type for the thermostat field
Field descriptions
Δt::Tthe time between timesteps, assumed to be in atomic unitssteps::Integerthe number of timesteps for the simulationt₀::Tthe starting time of the simulation, assumed to be in atomic units; defaults to 0thermostat::Cthe thermostat for the simulation; many options are defined byNBodySimulator, but a user could also define a custom thermostat; defaults to theNullThermostatsimulator::Sthe algorithm to be used for the ODE; defaults to VelocityVerlet
Atomistic.UnimplementedError — TypeUnimplementedErrorException thrown in default implementation of API to indicate that an implementator did not provide an implementation of a particular API function.
Atomistic.get_boundary_conditions — Methodget_boundary_conditions(result::MolecularDynamicsResult)::SVector{3, BoundaryCondition}Extract the boundary conditions of the simulation from the simulation result in an AtomBase-compatible format.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.get_bounding_box — Methodget_bounding_box(result::MolecularDynamicsResult)::SVector{3, SVector{3, Unitful.Length}}Extract the unit-anotated bounding box of the simulation from the simulation result in an AtomBase-compatible format.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.get_num_bodies — Methodget_num_bodies(result::MolecularDynamicsResult)::IntegerExtract the number of bodies in the simulation from the simulation result.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.get_particles — Methodget_particles(result::MolecularDynamicsResult, t::Integer)::AbstractVector{Atom}Extract the position of each particle in the system at a particular timestep from the simulation result in an AtomBase-compatible format. Automatically defaults to the end of the simulation when t is not passed.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.get_positions — Methodget_positions(result::MolecularDynamicsResult, t::Integer)::AbstractVector{SVector{3, Unitful.Length}}Extract the unit-anotated position of each particle in the system at a particular timestep from the simulation result in an AtomBase-compatible format. The returned positions should be normalized such that they lie within the bounding box of the system, even if the system is periodic. Automatically defaults to the end of the simulation when t is not passed.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.get_system — Methodget_system(result::MolecularDynamicsResult, t::Integer)::AbstractSystem{3}Extract the underlying system at a particular timestep from the simulation result. Automatically defaults to the end of the simulation when t is not passed.
The default implementation combines the results of get_particles, get_bounding_box, get_boundary_conditions, and get_time at the timestep to create a FlexibleSystem wrapped in a DyanmicSystem. An implementer of this API could implement a method of this function for their custom result type if it supports a more efficient way to extract the system.
Atomistic.get_time — Methodget_time(result::MolecularDynamicsResult, t::Integer)::Unitful.TimeExtract the unit-anotated time of the simulation at a particular timestep from the simulation result. Automatically defaults to the end of the simulation when t is not passed.
The default implementation extracts the time from the result of get_time_range.
Atomistic.get_time_range — Methodget_time_range(result::MolecularDynamicsResult)::AbstractVector{Unitful.Time}Extract the unit-anotated time range of the simulation from the simulation result.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.get_velocities — Methodget_velocities(result::MolecularDynamicsResult, t::Integer)::AbstractVector{SVector{3, Unitful.Velocity}}Extract the unit-anotated velocity of each particle in the system at a particular timestep from the simulation result in an AtomBase-compatible format. Automatically defaults to the end of the simulation when t is not passed.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.kinetic_energy — Methodkinetic_energy(result::MolecularDynamicsResult, t::Integer)::Unitful.EnergyExtract the unit-anotated kinetic energy of the simulation at a particular timestep from the simulation result. Automatically defaults to the end of the simulation when t is not passed.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.plot_energy! — Methodplot_energy!(p::Plot, result::MolecularDynamicsResult, stride::Integer; first_plot::Bool = false, time_unit = u"ps", energy_unit = u"hartree")::PlotPlot the kinetic, potential, and total energy of a MolecularDynamicsResult against time, sampling every stride points. Converts time and energy quantities to the specified units. Add the new lines to an existing plot and update the legend only if it is the first plot. If it is not the first plot, add a vertical line to differentiate the segments of the simulation in the plot.
Atomistic.plot_energy — Methodplot_energy(result::MolecularDynamicsResult, stride::Integer; time_unit = u"ps", energy_unit = u"hartree")::PlotPlot the kinetic, potential, and total energy of a MolecularDynamicsResult against time, sampling every stride points. Converts time and energy quantities to the specified units.
Atomistic.plot_rdf — Functionplot_rdf(result::MolecularDynamicsResult, σ::Unitful.Length, start::Integer = 1, stop::Integer = length(result))::Plot
Plot the radial distribution function of a `MolecularDynamicsResult` averaging over the timesteps in `start:stop`
Use `σ` (from Lennard Jones) as a normalization factor for the radius.Atomistic.plot_rdf — Functionplot_rdf(result::MolecularDynamicsResult, σ::Real, start::Integer = 1, stop::Integer = length(result))::PlotPlot the radial distribution function of a MolecularDynamicsResult averaging over the timesteps in start:stop Use σ (from Lennard Jones) as a normalization factor for the radius.
Atomistic.plot_temperature! — Methodplot_temperature!(p::Plot, result::MolecularDynamicsResult, stride::Integer; first_plot::Bool = false; first_plot::Bool = false, time_unit = u"ps", temperature_unit = u"K")::PlotPlot the temperature of a MolecularDynamicsResult against time, sampling every stride points. Converts time and temperature quantities to the specified units. Adds the new line to an existing plot and update the legend only if it is the first plot. If it is not the first plot, adds a vertical line to differentiate the segments of the simulation in the plot.
Atomistic.plot_temperature — Methodplot_temperature(result::MolecularDynamicsResult, stride::Integer; time_unit = u"ps", temperature_unit = u"K")::PlotPlot the temperature of a MolecularDynamicsResult against time, sampling every stride points. Converts time and temperature quantities to the specified units.
Atomistic.potential_energy — Methodpotential_energy(result::MolecularDynamicsResult, t::Integer)::Unitful.EnergyExtract the unit-anotated potential energy of the simulation at a particular timestep from the simulation result. Automatically defaults to the end of the simulation when t is not passed.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.rdf — Methodrdf(result::MolecularDynamicsResult, start::Integer, stop::Integer)::Tuple{AbstractVector{Unitful.Length},AbstractVector{Real}}Calculate the radial distribution function from the simulation result.
To include only a portion of the timesteps for faster computation, increase the start parameter or decrease the stop parameter; by default, the full time range will be iterated over. The result is a named tuple of vectors: - r: unit-annotated interparticle radial distance bins - g: distribution value of each bin
The default implementation provided uses the provided implementations of get_num_bodies, get_bounding_box, get_boundary_conditions, and get_positions. An implementor of the API could use a built in implementation if one that fits this spec is available.
Atomistic.reference_temperature — Methodreference_temperature(result::MolecularDynamicsResult)::Union{Unitful.Temperature,Missing}Extract the unit-anotated reference temperature of the simulation from the simulation result. If there is no thermostat with a reference temperature in this simulation, return missing.
An implementer of this API should implement a method of this function for their custom result type if it supports thermostats. If not implmented, the default implemenation just returns missing.
Atomistic.simulate — Methodsimulate(system::AbstractSystem{3}, simulator::MolecularDynamicsSimulator, potential::AbstractPotential)::MolecularDynamicsResultRun a molecular dynamics simulation configured with a particular simulator and potential with any abstract system.
An implementer of this API should implement a method of this function for their custom simulator type. If the simulator has a fast path for some types of potential, those should be implemented with multiple dispatch.
Atomistic.temperature — Methodtemperature(result::MolecularDynamicsResult, t::Integer)::Unitful.TemperatureExtract the unit-anotated temperature of the simulation at a particular timestep from the simulation result. Automatically defaults to the end of the simulation when t is not passed.
An implementer of this API should implement a method of this function for their custom result type.
Atomistic.total_energy — Methodtotal_energy(result::MolecularDynamicsResult, t::Integer)::Unitful.EnergyExtract the unit-anotated total energy of the simulation at a particular timestep from the simulation result. Automatically defaults to the end of the simulation when t is not passed.
The default implementation simply sums the results of kinetic_energy and potential_energy at the timestep. An implementer of this API could implement a method of this function for their custom result type if it supports a more efficient way to calculate this quantity.
Molly.visualize — Methodvisualize(result::MollyResult, filename::String; kwargs...)Animate a MollyResult using Molly's visualize function and store the result in a file. Any kwargs will be passed directly to the visualize function.
This function is only available when GLMakie is imported.
RecipesBase.animate — Methodanimate(result::NBSResult, filename::String; kwargs...)Animate an NBSResult using the interface provided by NBodySimulator.jl and store the result in a file.