Molly API
The API reference can be found here.
Molly also re-exports StaticArrays.jl and Unitful.jl, making the likes of SVector
and 1.0u"nm"
available when you call using Molly
.
The visualize
function is only available once you have called using Makie
. Requires.jl is used to lazily load this code when Makie.jl is available, which reduces the dependencies of the package.
Molly.AbstractCutoff
Molly.AndersenThermostat
Molly.Angletype
Molly.Atom
Molly.AtomMin
Molly.Atomtype
Molly.Bondtype
Molly.CoordinateLogger
Molly.Coulomb
Molly.DistanceNeighborFinder
Molly.EnergyLogger
Molly.GeneralInteraction
Molly.Gravity
Molly.HarmonicAngle
Molly.HarmonicBond
Molly.Interaction
Molly.LennardJones
Molly.Logger
Molly.Mie
Molly.NeighborFinder
Molly.NoCutoff
Molly.NoNeighborFinder
Molly.NoThermostat
Molly.ShiftedForceCutoff
Molly.ShiftedPotentialCutoff
Molly.Simulation
Molly.Simulator
Molly.SoftSphere
Molly.SpecificInteraction
Molly.StructureWriter
Molly.TemperatureLogger
Molly.Thermostat
Molly.Torsion
Molly.Torsiontype
Molly.TreeNeighborFinder
Molly.VelocityFreeVerlet
Molly.VelocityVerlet
Molly.accelerations
Molly.adjust_bounds
Molly.apply_thermostat!
Molly.displacements
Molly.distances
Molly.energy
Molly.find_neighbors!
Molly.force
Molly.log_property!
Molly.mass
Molly.maxwellboltzmann
Molly.placeatoms
Molly.rdf
Molly.readinputs
Molly.simulate!
Molly.temperature
Molly.vector
Molly.vector1D
Molly.velocity
Molly.visualize
Molly.AbstractCutoff
— TypeA general type of cutoff encoding the approximation used for a potential. Interactions can be parameterized by the cutoff behavior.
Molly.AndersenThermostat
— TypeAndersenThermostat(coupling_const)
Rescale random velocities according to the Andersen thermostat.
Molly.Angletype
— TypeAngletype(th0, cth)
Gromacs angle type.
Molly.Atom
— TypeAtom(; <keyword arguments>)
An atom and its associated information. Properties unused in the simulation or in analysis can be left with their default values. This type cannot be used on the GPU as it is not isbits
- use AtomMin
instead.
Arguments
attype::AbstractString=""
: the type of the atom.name::AbstractString=""
: the name of the atom.resnum::Integer=0
: the residue number if the atom is part of a polymer.resname::AbstractString=""
: the residue name if the atom is part of a polymer.charge::C=0.0u"q"
: the charge of the atom, used for electrostatic interactions.mass::M=0.0u"u"
: the mass of the atom.σ::S=0.0u"nm"
: the Lennard-Jones finite distance at which the inter-particle potential is zero.ϵ::E=0.0u"kJ * mol^-1"
: the Lennard-Jones depth of the potential well.
Molly.AtomMin
— TypeAtomMin(; <keyword arguments>)
An atom and minimal information required for simulation. Properties unused in the simulation or in analysis can be left with their default values. This type is isbits
and can be used on the GPU - use Atom
in other contexts to store more information.
Arguments
charge::C=0.0u"q"
: the charge of the atom, used for electrostatic interactions.mass::M=0.0u"u"
: the mass of the atom.σ::S=0.0u"nm"
: the Lennard-Jones finite distance at which the inter-particle potential is zero.ϵ::E=0.0u"kJ * mol^-1"
: the Lennard-Jones depth of the potential well.
Molly.Atomtype
— TypeAtomtype(charge, mass, σ, ϵ)
Gromacs atom type.
Molly.Bondtype
— TypeBondtype(b0, kb)
Gromacs bond type.
Molly.CoordinateLogger
— TypeCoordinateLogger(n_steps; dims=3)
Log the coordinates throughout a simulation.
Molly.Coulomb
— TypeCoulomb(; coulomb_const, cutoff, nl_only, force_unit, energy_unit)
The Coulomb electrostatic interaction.
Molly.DistanceNeighborFinder
— TypeDistanceNeighborFinder(nb_matrix, n_steps, dist_cutoff)
Find close atoms by distance.
Molly.EnergyLogger
— TypeEnergyLogger(n_steps)
Log the energy of the system throughout a simulation.
Molly.GeneralInteraction
— TypeA general interaction that will apply to all or most atom pairs. Custom general interactions should sub-type this type.
Molly.Gravity
— TypeGravity(; G, nl_only)
The gravitational interaction.
Molly.HarmonicAngle
— TypeHarmonicAngle(; i, j, k, th0, cth)
A bond angle between three atoms.
Molly.HarmonicBond
— TypeHarmonicBond(; i, j, b0, kb)
A harmonic bond between two atoms.
Molly.Interaction
— TypeAn interaction between atoms that contributes to forces on the atoms.
Molly.LennardJones
— TypeLennardJones(; cutoff, nl_only, force_unit, energy_unit, skip_shortcut)
The Lennard-Jones 6-12 interaction. The potential is given by
\[V_{ij}(r_{ij}) = 4\varepsilon_{ij} \left[\left(\frac{\sigma_{ij}}{r_{ij}}\right)^{12} - \left(\frac{\sigma_{ij}}{r_{ij}}\right)^{6}\right]\]
and the force on each atom by
\[\begin{aligned} \vec{F}_i &= 24\varepsilon_{ij} \left(2\frac{\sigma_{ij}^{12}}{r_{ij}^{13}} - \frac{\sigma_{ij}^6}{r_{ij}^{7}}\right) \frac{\vec{r}_{ij}}{r_{ij}} \\ &= \frac{24\varepsilon_{ij}}{r_{ij}^2} \left[2\left(\frac{\sigma_{ij}^{6}}{r_{ij}^{6}}\right)^2 -\left(\frac{\sigma_{ij}}{r_{ij}}\right)^{6}\right] \vec{r}_{ij} \end{aligned}\]
Molly.Logger
— TypeA way to record a property, e.g. the temperature, throughout a simulation. Custom loggers should sub-type this type.
Molly.Mie
— TypeMie(; m, n, cutoff, nl_only, force_unit, energy_unit, skip_shortcut)
The Mie generalized interaction. When m
equals 6 and n
equals 12 this is equivalent to the Lennard-Jones interaction.
Molly.NeighborFinder
— TypeA way to find near atoms to save on simulation time. Custom neighbor finders should sub-type this type.
Molly.NoCutoff
— TypeNoCutoff()
Placeholder cutoff that does not alter forces or potentials.
Molly.NoNeighborFinder
— TypeNoNeighborFinder()
Placeholder neighbor finder that returns no neighbors. When using this neighbor finder, ensure that nl_only
for the interactions is set to false
.
Molly.NoThermostat
— TypeNoThermostat()
Placeholder thermostat that does nothing.
Molly.ShiftedForceCutoff
— TypeShiftedForceCutoff(cutoff_dist)
Cutoff that shifts the force to be continuous at a specified cutoff point.
Molly.ShiftedPotentialCutoff
— TypeShiftedPotentialCutoff(cutoff_dist)
Cutoff that shifts the potential to be continuous at a specified cutoff point.
Molly.Simulation
— TypeSimulation(; <keyword arguments>)
The data needed to define and run a molecular simulation. Properties unused in the simulation or in analysis can be left with their default values.
Arguments
simulator::Simulator
: the type of simulation to run.atoms::A
: the atoms, or atom equivalents, in the simulation. Can be of any type.specific_inter_lists::SI=()
: the specific interactions in the simulation, i.e. interactions between specific atoms such as bonds or angles. Typically aTuple
.general_inters::GI=()
: the general interactions in the simulation, i.e. interactions between all or most atoms such as electrostatics. Typically aTuple
.coords::C
: the coordinates of the atoms in the simulation. Typically aVector
ofSVector
s of any dimension and typeT
, whereT
is anAbstractFloat
type.velocities::V=zero(coords)
: the velocities of the atoms in the simulation, which should be the same type as the coordinates. The meaning of the velocities depends on the simulator used, e.g. for theVelocityFreeVerlet
simulator they represent the previous step coordinates for the first step.temperature::T=0.0
: the temperature of the simulation.box_size::B
: the size of the cube in which the simulation takes place.neighbor_finder::NeighborFinder=NoNeighborFinder()
: the neighbor finder used to find close atoms and save on computation.thermostat::Thermostat=NoThermostat()
: the thermostat which applies during the simulation.loggers::Dict{String, <:Logger}=Dict()
: the loggers that record properties of interest during the simulation.timestep::S=0.0
: the timestep of the simulation.n_steps::Integer=0
: the number of steps in the simulation.force_unit::F=u"kJ * mol^-1 * nm^-1"
: the unit of force used in the simulation.energy_unit::E=u"kJ * mol^-1"
: the unit of energy used in the simulation.gpu_diff_safe::Bool
: whether to use the GPU implementation. Defaults toisa(coords, CuArray)
.
Molly.Simulator
— TypeA type of simulation to run, e.g. leap-frog integration or energy minimisation. Custom simulators should sub-type this type.
Molly.SoftSphere
— TypeSoftSphere(; cutoff, nl_only, force_unit, energy_unit, skip_shortcut)
The soft-sphere potential.
Molly.SpecificInteraction
— TypeA specific interaction between sets of specific atoms, e.g. a bond angle. Custom specific interactions should sub-type this type.
Molly.StructureWriter
— TypeStructureWriter(n_steps, filepath)
Write 3D output structures to the PDB file format throughout a simulation.
Molly.TemperatureLogger
— TypeTemperatureLogger(n_steps)
TemperatureLogger(T, n_steps)
Log the temperature throughout a simulation.
Molly.Thermostat
— TypeA way to keep the temperature of a simulation constant. Custom thermostats should sub-type this type.
Molly.Torsion
— TypeTorsion(; i, j, k, l, f1, f2, f3, f4)
A dihedral torsion angle between four atoms.
Molly.Torsiontype
— TypeTorsiontype(f1, f2, f3, f4)
Gromacs torsion type.
Molly.TreeNeighborFinder
— TypeTreeNeighborFinder(nb_matrix, n_steps, dist_cutoff)
Find close atoms by distance using a tree search.
Molly.VelocityFreeVerlet
— TypeVelocityFreeVerlet()
The velocity-free Verlet integrator, also known as the Störmer method. In this case the velocities
given to the Simulator
act as the previous step coordinates for the first step.
Molly.VelocityVerlet
— TypeVelocityVerlet()
The velocity Verlet integrator.
Molly.accelerations
— Methodaccelerations(simulation; parallel=true)
Calculate the accelerations of all atoms using the general and specific interactions and Newton's second law.
Molly.adjust_bounds
— Methodadjust_bounds(c, box_size)
Ensure a coordinate is within the simulation box and return the coordinate.
Molly.apply_thermostat!
— Methodapply_thermostat!(simulation, thermostat)
Apply a thermostat to modify a simulation. Custom thermostats should implement this function.
Molly.displacements
— Methoddisplacements(coords, box_size)
Get the pairwise vector displacements of a set of coordinates, accounting for the periodic boundary conditions.
Molly.distances
— Methoddistances(coords, box_size)
Get the pairwise distances of a set of coordinates, accounting for the periodic boundary conditions.
Molly.energy
— Methodenergy(s)
Compute the total energy of the system.
Molly.find_neighbors!
— Methodfind_neighbors!(simulation, neighbor_finder, step_n; parallel=true)
Obtain a list of close atoms in a system. Custom neighbor finders should implement this function.
Molly.force
— Functionforce(inter, coord_i, coord_j, atom_i, atom_j, box_size)
Calculate the force between a pair of atoms due to a given interation type. Custom interaction types should implement this function.
Molly.log_property!
— Methodlog_property!(logger, simulation, step_n)
Log a property thoughout a simulation. Custom loggers should implement this function.
Molly.mass
— Methodmass(atom)
The mass of an atom.
Molly.maxwellboltzmann
— Methodmaxwellboltzmann(mass, temperature)
Draw a speed along one dimension in accordance with the Maxwell-Boltzmann distribution.
Molly.placeatoms
— Methodplaceatoms(n_atoms, box_size, min_dist; dims=3)
Obtain n_atoms
3D coordinates in a cube of length box_size
where no two points are closer than min_dist
, accounting for periodic boundary conditions.
Molly.rdf
— Methodrdf(coords, box_size; npoints=200)
Get the radial distribution function of a set of coordinates. This describes how density varies as a function of distance from each atom. Returns a list of distance bin centres and a list of the corresponding densities.
Molly.readinputs
— Methodreadinputs(topology_file, coordinate_file; atom_min=false, units=true)
readinputs(T, topology_file, coordinate_file; atom_min=false, units=true)
Read a Gromacs topology flat file, i.e. all includes collapsed into one file, and a Gromacs coordinate file. Returns the atoms, specific interaction lists, general interaction lists, non-bonded matrix, coordinates and box size. atom_min
determines whether the returned atoms consist of the GPU-compatible AtomMin
or the more verbose but GPU-incompatible Atom
. units
determines whether the returned values have units.
Molly.simulate!
— Methodsimulate!(simulation; parallel=true)
simulate!(simulation, n_steps; parallel=true)
simulate!(simulation, simulator, n_steps; parallel=true)
Run a simulation according to the rules of the given simulator. Custom simulators should implement this function.
Molly.temperature
— Methodtemperature(simulation)
Calculate the temperature of a system from the kinetic energy of the atoms.
Molly.vector
— Methodvector(c1, c2, box_size)
Displacement between two coordinate values, accounting for the bounding box. The minimum image convention is used, so the displacement is to the closest version of the coordinates accounting for the periodic boundaries.
Molly.vector1D
— Methodvector1D(c1, c2, box_size)
Displacement between two 1D coordinate values, accounting for the bounding box. The minimum image convention is used, so the displacement is to the closest version of the coordinate accounting for the periodic boundaries.
Molly.velocity
— Methodvelocity(mass, temperature; dims=3)
Generate a random velocity from the Maxwell-Boltzmann distribution.
Molly.visualize
— Functionvisualize(coord_logger, box_size, out_filepath; <keyword arguments>)
Visualize a simulation as an animation. This function is only available when Makie is imported. It can take a while to run, depending on the length and size of the simulation.
Arguments
connections=Tuple{Int, Int}[]
: pairs of atoms indices to link with bonds.connection_frames
: the frames in which bonds are shown. Is a list of the same length as the number of frames, where each item is a list ofBool
s of the same length asconnections
. Defaults to alwaystrue
.trails::Integer=0
: the number of preceding frames to show as transparent trails.framerate::Integer=30
: the frame rate of the animation.color=:purple
: the color of the atoms. Can be a single color or a list of colors of the same length as the number of atoms.connection_color=:orange
: the color of the bonds. Can be a single color or a list of colors of the same length asconnections
.markersize=0.1
: the size of the atom markers.linewidth=2.0
: the width of the bond lines.transparency=true
: whether transparency is active on the plot.kwargs...
: other keyword arguments are passed to the plotting function.