Getting started
To get started with the package we can identify the necessary ingredients to perform a simple classical dynamics simulation and walk through how to set up the simulation.
Atoms
First, we must define the particles in the simulation. For this purpose we provide the Atoms
type which will contain the symbols, atomic numbers and masses for our atoms. Technically these need not be actual atoms and be a generic particle.
If using real atoms, then they can be constructed using the chemical symbols as a Vector
of Julia's Symbol
types, a Vector{Symbol}
:
julia> Atoms([:H, :C])
Atoms{Float64}([:H, :C], [1, 6], [1837.4715941070515, 21894.713607956142])
You can see that this contains two atoms labelled by their atomic numbers with their masses in atomic units.
Internally atomic units are used for all quantities. This makes things simple when performing nonadiabatic dynamics. Unitful.jl and UnitfulAtomic.jl can be used to help with unit transformations, and many functions will directly accept Unitful quantities and handle the conversions for you.
Alternatively, if not using real atoms, Atoms
can be created using a Vector{<:Real}
where the provided numbers are the masses of the particles.
julia> Atoms([1, 2, 3, 4, 5, 6])
Atoms{Float64}([:X, :X, :X, :X, :X, :X], [0, 0, 0, 0, 0, 0], [1.0, 2.0, 3.0, 4.0, 5.0, 6.0])
A more detailed look into the Atoms
type along with a description of how to save and load structures can be found here.
Representing atomic positions and velocities
This package chooses to separate the dynamical variables from the static atomic parameters included in the Atoms
type. This allows us to easily interface with other numerical packages like DifferentialEquations.jl and AdvancedMH.jl. As such, both positions and velocities are represented using Julia's standard Array
type, specifically as an Array{T,2}
or the Matrix{T}
type, which are equivalent. If you are new to Julia, you can find a description of the Array
here. The first dimension contains each atomic degree of freedom, and the second dimension contains each atom. For example, a 3D system with two atoms would have positions:
julia> using Symbolics
julia> @variables x1, y1, z1, x2, y2, z2
6-element Vector{Symbolics.Num}: x1 y1 z1 x2 y2 z2
julia> r = [x1 x2; y1 y2; z1 z2]
3×2 Matrix{Symbolics.Num}: x1 x2 y1 y2 z1 z2
Symbolics is a package available from the General registry. You will have to add it to your current environment using pkg> add Symbolics
to be able to reproduce this example. Throughout the documentation we occasionally use external packages, if you run into an error you will likely have to add
the package before being able to use it. Refer to the Julia manual for further information on installing packages.
For a 1D system it would be necessary to create a 1x1 matrix:
julia> r = fill(x1, (1,1))
1×1 Matrix{Symbolics.Num}: x1
Velocities are handled in the same way as positions and the data structures are the same. Usually manual initialisation like this will only be necessary for small model systems, whereas full dimensional model system will be read from a file instead. This is explored in the Atoms
documentation.
We can also perform simulations using ring polymers which have multiple replicas of each atom, these are implemented using Array{T,3}
where the third dimension is used for each ring polymer bead. For more information, see the ring polymer methods in the dynamics methods section.
Models
The next ingredient required to set up the simulation is the Model
, i.e., the potentials in which the system evolves. These Model
s are provided by NQCModels.jl, which is a convenient infrastructure for defining different kinds of models for adiabatic and nonadiabatic dynamics. These models can range from simple analytic potentials all the way up to multi-dimensional ab initio potentials. Refer to the NQCModels.jl page for information on the available models and a description of how to implement further models.
For now we can look at an AdiabaticModel
which provides a simple harmonic potential energy function.
julia> model = Harmonic()
Harmonic{Float64, Float64, Float64} m: Float64 1.0 ω: Float64 1.0 r₀: Float64 0.0 dofs: Int64 1
Here, the four parameters (m
, ω
, r₀
and dofs
) for this model are shown along with their types and default values. These values can be modified by specifying a new value in the constructor. For example for m
:
julia> model = Harmonic(m=0.4)
Harmonic{Float64, Float64, Float64} m: Float64 0.4 ω: Float64 1.0 r₀: Float64 0.0 dofs: Int64 1
Many of the models use Parameters.jl to provide convenient keyword constructors and formatted printing for the models. The Harmonic
model above is defined using the @with_kw
macro from Parameters.jl
to give it a set of default parameters. Each of these can be modified by specifying a new value using keyword arguments in the constructor as demonstrated above.
Adiabatic models implement two functions to calculate the total energy and the forces, respectively: potential(model, R)
and derivative(model, R)
.
Let's try these out and take a look at the results:
julia> potential(model, hcat(25.0))
125.0
julia> derivative(model, hcat(25.0))
1×1 Matrix{Float64}: 10.0
All models accept an R::AbstractMatrix
for the argument representing the positions of the particles in the system. These are structured such that size(R) = (dofs, natoms)
where dofs
is the number of degrees of freedom for each atom, and natoms
is the number of atoms in the simulation.
Since this is a 1D model, we use hcat
to quickly create a 1x1 matrix.
To make sure the model is what we expect, we can plot the potential and derivative using a custom plotting recipe. This looks pretty harmonic to me!
using Plots
plot(-5:0.1:5, model)
Plotting recipes currently only exist for 1D models. For more complex models you will have to handle the plotting manually.
Simulation
To control all simulation parameters in one environment, we use the Simulation
type which will contain both the Atoms
and Model
s explained above, along with any extra information required for the simulation.
julia> sim = Simulation{Classical}(Atoms(:H), model)
Simulation{Classical}: Atoms{Float64}([:H], [1], [1837.4715941070515]) Harmonic{Float64, Float64, Float64} m: Float64 0.4 ω: Float64 1.0 r₀: Float64 0.0 dofs: Int64 1
Here, we have specified that each atom has a single degree of freedom and have not provided a simulation cell. Classical
is a type parameter, and specifies the dynamics method that we want to use. Check out Dynamics methods to learn about the other kinds of dynamics available.
Technically Simulation(atoms, model)
is equivalent to Simulation{Classical}(atoms, model)
since Classical
is the default.
Dynamics variables
The final ingredient before we can perform our simulation is the initial positions and velocities of our particles. For each dynamics type, the method DynamicsVariables
is implemented and creates the dynamics variables for us. For classical dynamics we must provide a Matrix
of velocities and of positions. These should have size = (dofs, natoms)
, matching the arguments of the potential
and derivative
functions. Usually the initial coordinates would have some physical significance, perhaps sampled from a relevant distribution, but here we use random numbers for simplicity.
julia> v = rand(3, 3);
julia> r = rand(3, 3);
julia> DynamicsVariables(sim, v, r)
ComponentVector{Float64}(v = [0.16225348864938283 0.9829777922862497 0.5912386603731087; 0.047619934910227424 0.214638250946262 0.9953523072299958; 0.41103666443749043 0.8742949983221534 0.6675813983315546], r = [0.16660924665147503 0.9857149711909038 0.6909355873337776; 0.8087109347835648 0.47065352599937094 0.6891415809298602; 0.6024746390763626 0.12319794038862919 0.3440136317275013])
Since DifferentialEquations.jl
requires AbstractArray
s for the dynamics variables, we use ComponentArrays.jl which allow us to conveniently store all the required information for different types of dynamics.
Bringing it all together
We have now covered all the parts necessary to perform our first classical dynamics simulation. Let's quickly set up our simulation parameters using what we've learned. Here we'll have two atoms in a harmonic potential, each with a single degree of freedom.
julia> atoms = Atoms([:H, :C])
Atoms{Float64}([:H, :C], [1, 6], [1837.4715941070515, 21894.713607956142])
julia> sim = Simulation{Classical}(atoms, Harmonic(ω=50.0))
Simulation{Classical}: Atoms{Float64}([:H, :C], [1, 6], [1837.4715941070515, 21894.713607956142]) Harmonic{Float64, Float64, Float64} m: Float64 1.0 ω: Float64 50.0 r₀: Float64 0.0 dofs: Int64 1
julia> z = DynamicsVariables(sim, randn(size(sim)), randn(size(sim)))
ComponentVector{Float64}(v = [-0.8656368220680416 -1.5703832464742633], r = [-1.6853073306013029 0.3151755924118576])
Now, we can finally run the trajectory using the run_dynamics
function. This takes three positional arguments: the simulation parameters sim, the time span we want to solve for
tspan, and the dynamics variables
z. For classical dynamics we also provide a timestep
dtsince we're using the
VelocityVerlet` algorithm by default.
Each method will default to an appropriate integration algorithm though it is possible to specify via a keyword argument to run_dynamics
if an alternative algorithm is preferred. Refer to the dynamics documentation for more information.
The final keyword argument output
is used to specify the quantities we want to save during the dynamics. A list of the available quantities can be found here.
run_dynamics
returns a Dictionary
from Dictionaries.jl
that has entries containing the time and the output quantities saved at each time step.
tspan = (0.0, 50.0)
solution = run_dynamics(sim, (0.0, 50.0), z;
dt=0.1, output=(OutputPosition, OutputVelocity))
3-element Dictionaries.Dictionary{Symbol, Any}
:Time │ [0.0, 0.1, 0.2, 0.30000000000000004, 0.4, 0.5, 0.6, 0.7, 0.7…
:OutputPosition │ [[-1.6853073306013029 0.3151755924118576], [-1.7604061598546…
:OutputVelocity │ [[-0.8656368220680416 -1.5703832464742633], [-0.631230920127…
Here you can see the output containing the time steps and the output quantities we specified. These can be accessed directly as shown here:
julia> solution[:Time]
501-element Vector{Float64}: 0.0 0.1 0.2 0.30000000000000004 0.4 0.5 0.6 0.7 0.7999999999999999 0.8999999999999999 ⋮ 49.20000000000043 49.30000000000043 49.40000000000043 49.50000000000043 49.600000000000435 49.700000000000436 49.80000000000044 49.90000000000044 50.0
julia> solution[:OutputPosition]
501-element Vector{Matrix{Float64}}: [-1.6853073306013029 0.3151755924118576] [-1.760406159854691 0.15795732958502906] [-1.8115535146267798 0.0005587066351893358] [-1.8380535017819133 -0.15684055426160554] [-1.839545571703432 -0.3140607302008737] [-1.8160094238064677 -0.4709223027623764] [-1.7677652827406487 -0.6272461629893006] [-1.6954695415247916 -0.7828538158999057] [-1.6001058308915794 -0.9375675842981146] [-1.482971636349619 -1.0912108116503363] ⋮ [-1.6529746265773047 3.5008062513018188] [-1.547089295385074 3.6026564353746635] [-1.4201548048432953 3.7003930049697358] [-1.273898181521371 3.793904361722077] [-1.1103093422129944 3.883083731733629] [-0.9316140198197382 3.967829287490672] [-0.7402434808026237 4.048044264133348] [-0.5388014462155005 4.123637069944496] [-0.33002866638209816 4.194521390931357]
As with the models, we provide custom plotting recipes to quickly visualise the results before performing further analysis by manually accessing the fields of the solution table. To use these recipes, simply provide the solution to the plot
function from Plots.jl
and give the name of the output quantity as the second argument. This will only work if this quantity was specified in run_dynamics
.
plot(solution, :OutputPosition)
plot!(solution, :OutputVelocity)
Ensemble simulations
We have shown how to perform a single trajectory, but usually we are interested in performing many and calculating observables using statistical methods. Running more trajectories is as simple as providing the trajectories
keyword to run_dynamics
, but we'll go through this in more detail in the Ensemble simulations section.
What's next?
Now that we've covered the basics of classical dynamics, we're ready to explore the world of nonadiabatic dynamics. All the dynamics methods follow these patterns and anything you find elsewhere in the documentation should now seem relatively familiar.