API

Molly API

The API reference can be found here.

Molly also re-exports StaticArrays.jl, making the likes of SVector 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.

AndersenThermostat(coupling_const)

Rescale random velocities according to the Andersen thermostat.

Molly.AngletypeType.
Angletype(th0, cth)

Gromacs angle type.

Molly.AtomType.
Atom(; <keyword arguments>)

An atom and its associated information. Properties unused in the simulation or in analysis can be left with their default values.

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::T=0.0: the charge of the atom, used for electrostatic interactions.
  • mass::T=0.0: the mass of the atom.
  • σ::T=0.0: the Lennard-Jones finite distance at which the inter-particle potential is zero.
  • ϵ::T=0.0: the Lennard-Jones depth of the potential well.
Molly.AtomtypeType.
Atomtype(mass, charge, σ, ϵ)

Gromacs atom type.

Molly.BondtypeType.
Bondtype(b0, kb)

Gromacs bond type.

CoordinateLogger(n_steps; dims=3)

Log the coordinates throughout a simulation.

Molly.CoulombType.
Coulomb(nl_only)

The Coulomb electrostatic interaction.

DistanceNeighbourFinder(nb_matrix, n_steps, dist_cutoff)
DistanceNeighbourFinder(nb_matrix, n_steps)

Find close atoms by distance.

A general interaction that will apply to all atom pairs. Custom general interactions should sub-type this type.

Molly.GravityType.
Gravity(nl_only, G)

The gravitational interaction.

HarmonicAngle(i, j, k, th0, cth)

A bond angle between three atoms.

HarmonicBond(i, j, b0, kb)

A harmonic bond between two atoms.

An interaction between atoms that contributes to forces on the atoms.

LennardJones(nl_only)

The Lennard-Jones 6-12 interaction.

Molly.LoggerType.

A way to record a property, e.g. the temperature, throughout a simulation. Custom loggers should sub-type this type.

A way to find near atoms to save on simulation time. Custom neighbour finders should sub-type this type.

NoNeighbourFinder()

Placeholder neighbour finder that returns no neighbours. When using this neighbour finder, ensure that nl_only for the interactions is set to false.

NoThermostat()

Placeholder thermostat that does nothing.

Simulation(; <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::Vector{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 a Tuple.
  • general_inters::GI=(): the general interactions in the simulation, i.e. interactions between all or most atoms such as electrostatics. Typically a Tuple.
  • coords::C: the coordinates of the atoms in the simulation. Typically a Vector of SVectors of any dimension and type T, where T is an AbstractFloat type.
  • velocities::C=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 the VelocityFreeVerlet simulator they represent the previous step coordinates for the first step.
  • temperature::T=0.0: the temperature of the simulation.
  • box_size::T: the size of the cube in which the simulation takes place.
  • neighbour_finder::NeighbourFinder=NoNeighbourFinder(): the neighbour 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::T=0.0: the timestep of the simulation.
  • n_steps::Integer=0: the number of steps in the simulation.
  • n_steps_made::Vector{Int}=[]: the number of steps already made during the simulation. This is a Vector to allow the struct to be immutable.
Molly.SimulatorType.

A type of simulation to run, e.g. leap-frog integration or energy minimisation. Custom simulators should sub-type this type.

SoftSphere(nl_only)

The soft-sphere potential.

A specific interaction between sets of specific atoms, e.g. a bond angle. Custom specific interactions should sub-type this type.

StructureWriter(n_steps, filepath)

Write 3D output structures to the PDB file format throughout a simulation.

TemperatureLogger(n_steps)
TemperatureLogger(T, n_steps)

Log the temperature throughout a simulation.

A way to keep the temperature of a simulation constant. Custom thermostats should sub-type this type.

Molly.TorsionType.
Torsion(i, j, k, l, f1, f2, f3, f4)

A dihedral torsion angle between four atoms.

Torsiontype(f1, f2, f3, f4)

Gromacs torsion type.

VelocityFreeVerlet()

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.

VelocityVerlet()

The velocity Verlet integrator.

accelerations(simulation, neighbours; parallel=true)

Calculate the accelerations of all atoms using the general and specific interactions and Newton's second law.

adjust_bounds(c, box_size)

Ensure a coordinate is within the simulation box and return the coordinate.

apply_thermostat!(simulation, thermostat)

Apply a thermostat to modify a simulation. Custom thermostats should implement this function.

displacements(coords, box_size)

Get the pairwise vector displacements of a set of coordinates, accounting for the periodic boundary conditions.

Molly.distancesMethod.
distances(coords, box_size)

Get the pairwise distances of a set of coordinates, accounting for the periodic boundary conditions.

find_neighbours(simulation, current_neighbours, neighbour_finder, step_n; parallel=true)

Obtain a list of close atoms in a system. Custom neighbour finders should implement this function.

Molly.force!Function.
force!(forces, interaction, simulation, atom_i, atom_j)

Update the force for an atom pair in response to a given interation type. Custom interaction types should implement this function.

log_property!(logger, simulation, step_n)

Log a property thoughout a simulation. Custom loggers should implement this function.

maxwellboltzmann(mass, temperature)
maxwellboltzmann(T, mass, temperature)

Draw from the Maxwell-Boltzmann distribution.

Molly.rdfMethod.
rdf(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.readinputsMethod.
readinputs(topology_file, coordinate_file)
readinputs(T, topology_file, coordinate_file)

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.

Molly.simulate!Method.
simulate!(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.temperatureMethod.
temperature(simulation)

Calculate the temperature of a system from the kinetic energy of the atoms.

Molly.vectorMethod.
vector(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.vector1DMethod.
vector1D(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.velocityMethod.
velocity(mass, temperature; dims=3)
velocity(T, mass, temperature; dims=3)

Generate a random velocity from the Maxwell-Boltzmann distribution.

Molly.visualizeFunction.
visualize(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 of Bools of the same length as connections. Defaults to always true.
  • 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 as connections.
  • 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.