API

LAMMPS.FixExternalType
FixExternal(callback, lmp::LMP, name::String, group::String, ncall::Int, napply::Int)

Creates a fix in lammps that calls a Julia function callback every ncall during the simulation.

lammps commands

The following command is executed in LAMMPS when FixExternal is called in order to setup the fix:

fix <name> <group> external pf/callback <ncall> <napply>

The FixExternal object gets passed to the callback function, it contains the parameters passed to FixExternal as fields:

  • lmp::LMP: The LAMMPS object.
  • name::String: The name of the fix.
  • group::String: The group of atoms to which the fix is applied.
  • ncall::Int: The number of timesteps between calls to the callback function.
  • napply::Int: The number of times the callback function is applied to the atoms in the group.

Additionally, the following fields are set by the callback function:

  • timestep::Int: The current timestep of the simulation.
  • nlocal::Int: The number of local atoms on the MPI rank.
  • ids::UnsafeArray{Int32, 1}: The IDs of the local atoms in an nall-element Vector.
  • x::UnsafeArray{Float64, 2}: The positions of the local atoms in a (3 × nall) Matrix.
  • f::UnsafeArray{Float64, 2}: The forces on the local atoms in a (3 × nlocal) Matrix.

These values are only valid during the callback execution and should not be used outside of it, as they refer directly to LAMMPS internal memory.

Here, the intention is for the callback to write forces to the f field, which will be applied to the atoms in the group every napply timesteps. f is not zeroed before the callback is called, so the forces from previous calls are preserved.

ghost atoms

FixExternal is invoked at the post_force stage of each timestep. Modified data of the ghost atoms will not be communicated back to the processor owning the respective atom, as at this point the reverse_comm stage of the timestep has already happened.

Contributions to the per-atom energies or virials can be set with set_energy! and set_virial!, respectively.

source
LAMMPS.InteractionConfigType
InteractionConfig(;
    system::Type{<:NamedTuple} = @NamedTuple{},
    atom::Type{<:NamedTuple} = @NamedTuple{type::Int32},
    backend::Union{Nothing, AbstractADType} = nothing
)

A configuration struct that encapsulates the system properties, atom properties, and an optional differentiation backend for use in external pair interactions.

Fields

  • system: A NamedTuple type that defines the global system properties to be extracted from LAMMPS. Defaults to an empty NamedTuple. A list of the available global properties can be found in the LAMMPS documentation
  • atom: A NamedTuple type that defines the per-atom properties to be extracted from LAMMPS. Defaults to @NamedTuple{type::Int32}. A list of the available atom properties can be found in the LAMMPS documentation
  • backend: An optional automatic differentiation backend used to automatically calculate interaction forces from a potential function. Defaults to nothing. The backends are defined through the ADTypes package.

Example

config = InteractionConfig(
    system = @NamedTuple{qqrd2e::Float64, boltz::Float64},
    atom = @NamedTuple{type::Int32, q::Float64},
    backend = nothing
)
source
LAMMPS.LMPType
LMP(f::Function, args=String[], comm=MPI.COMM_WORLD)

Create a new LAMMPS instance and call f on that instance while returning the result from f.

source
LAMMPS.LMPType
LMP(args::Vector{String}=String[], comm::MPI.Comm=MPI.COMM_WORLD)

Create a new LAMMPS instance while passing in a list of strings as if they were command-line arguments for the LAMMPS executable.

A full ist of command-line options can be found in the lammps documentation.

source
LAMMPS.PairExternalMethod
PairExternal(compute_potential, lmp::LMP, config::InteractionConfig, cutoff::Float64)

Defines a custom pair style in LAMMPS using a user-provided potential function.

limitations

PairExternal is implemented as a FixExternal and only works with newton off. Fixes that operate on forces, energies, or virials should also be defined after calling PairExternal. Furthermore, it is not possible to define multiple PairExternals for a simulation. Attempting to do this will overwrite the previously defined PairExternal.

lammps commands

The following commands are executed in LAMMPS during the setup process of PairExternal:

fix pair_julia all external pf/callback 1 1
pair_style zero <cutoff> nocoeff
pair_coeff * *
newton off <newton_bond> # turns off newton_pair while leaving newton_bond unchanged

compute_potential should have the following signature:

function compute_potential(r::Real, system::config.system, iatom::config.atom, jatom::config.atom)

If no differentiation backend is provided through config, the function should return a tuple of the pair energy and force magnitude. Otherwise it should return only the energy.

Examples

# without automatic differentiation
config = InteractionConfig(
    system = @NamedTuple{qqrd2e::Float64},
    atom = @NamedTuple{q::Float64},
)

PairExternal(lmp, config, 2.5) do r, system, iatom, jatom
    energy = system.qqrd2e * (iatom.q * jatom.q) / r
    force = system.qqrd2e * (iatom.q * jatom.q) / r^2
    return energy, force
end

# with automatic differentiation
import ADTypes: AutoEnzyme()
import Enzyme

config = InteractionConfig(
    system = @NamedTuple{qqrd2e::Float64},
    atom = @NamedTuple{q::Float64},
    backend = AutoEnzyme(),
)

PairExternal(lmp, config, 2.5) do r, system, iatom, jatom
    system.qqrd2e * (iatom.q * jatom.q) / r
end
source
LAMMPS.commandMethod
command(lmp::LMP, cmd::Union{String, Array{String}})

Process LAMMPS input commands from a String or from an Array of Strings.

A full list of commands can be found in the lammps documentation.

This function processes a multi-line string similar to a block of commands from a file. The string may have multiple lines (separated by newline characters) and also single commands may be distributed over multiple lines with continuation characters (’&’). Those lines are combined by removing the ‘&’ and the following newline character. After this processing the string is handed to LAMMPS for parsing and executing.

Arrays of Strings get concatenated into a single String inserting newline characters as needed.

LAMMPS.jl 0.4.1

Multiline string support """ and support for array of strings was added. Prior versions of LAMMPS.jl ignore newline characters.

Examples

LMP(["-screen", "none"]) do lmp
    command(lmp, """
        atom_modify map yes
        region cell block 0 2 0 2 0 2
        create_box 1 cell
        lattice sc 1
        create_atoms 1 region cell
        mass 1 1

        group a id 1 2 3 5 8
        group even id 2 4 6 8
        group odd id 1 3 5 7
    """)
end
source
LAMMPS.compute_neighborlistMethod
compute_neighborlist(lmp::LMP, id::String; request=0)

Retrieve neighbor list requested by a compute.

The neighbor list request from a compute is identified by the compute ID and the request ID. The request ID is typically 0, but will be > 0 in case a compute has multiple neighbor list requests.

Each neighbor list contains vectors of local indices of neighboring atoms. These can be used to index into Arrays returned from extract_atom.

source
LAMMPS.create_atomsMethod
create_atoms(
    lmp::LMP, x::Matrix{Float64}, id::Vector{Int32}, types::Vector{Int32};
    v::Union{Nothing,Matrix{Float64}}=nothing,
    image::Union{Nothing,Vector{IMAGEINT}}=nothing,
    bexpand::Bool=false
)

Create atoms for a LAMMPS instance. x contains the atom positions and should be a 3 by n Matrix{Float64}, where n is the number of atoms. id contains the id of each atom and should be a length n Vector{Int32}. types contains the atomic type (LAMMPS number) of each atom and should be a length n Vector{Int32}. v contains the associated velocities and should be a 3 by n Matrix{Float64}. image contains the image flags for each atom and should be a length n Vector{IMAGEINT}. bexpand is a Bool that defines whether or not the box should be expanded to fit the input atoms (default not).

source
LAMMPS.extract_atomMethod
extract_atom(lmp::LMP, name::String, lmp_type::_LMP_DATATYPE; with_ghosts=false)

Extract per-atom data from the lammps instance.

valid values for lmp_type:resulting return type:
LAMMPS_INTUnsafeArray{Int32, 1}
LAMMPS_INT_2DUnsafeArray{Int32, 2}
LAMMPS_DOUBLEUnsafeArray{Float64, 1}
LAMMPS_DOUBLE_2DUnsafeArray{Float64, 2}
LAMMPS_INT64UnsafeArray{Int64, 1}
LAMMPS_INT64_2DUnsafeArray{Int64, 2}
Info

The returned data may become invalid if a re-neighboring operation is triggered at any point after calling this method. If this has happened, trying to read from this data will likely cause julia to crash. To prevent this, copy the returned data

A table with suported name keywords can be found in the lammps documentation.

Arguments

  • with_ghosts: Determines wheter entries for ghost atoms are included. This is ignored for "mass", or when there is no ghost atom data available.
source
LAMMPS.extract_computeMethod
extract_compute(lmp::LMP, name::String, style::_LMP_STYLE_CONST, lmp_type::_LMP_TYPE)

Extract data provided by a compute command identified by the compute-ID. Computes may provide global, per-atom, or local data, and those may be a scalar, a vector or an array. Since computes may provide multiple kinds of data, it is required to set style and type flags representing what specific data is desired.

valid values for style:
STYLE_GLOBAL
STYLE_ATOM
STYLE_LOCAL
valid values for lmp_type:resulting return type:
TYPE_SCALARUnsafeArray{Float64, 0}
TYPE_VECTORUnsafeArray{Float64, 1}
TYPE_ARRAYUnsafeArray{Float64, 2}

Scalar values get returned as arrays with a single element. This way it's possible to modify the internal state of the LAMMPS instance even if the data is scalar.

Info

The returned data may become invalid as soon as another LAMMPS command has been issued at any point after calling this method. If this has happened, trying to read from this data will likely cause julia to crash. To prevent this, copy the returned data.

Examples

LMP(["-screen", "none"]) do lmp
    extract_compute(lmp, "thermo_temp", LMP_STYLE_GLOBAL, TYPE_SCALAR) |> println # [0.0]
end
source
LAMMPS.extract_globalMethod
extract_global(lmp::LMP, name::String, lmp_type::_LMP_DATATYPE)

Extract a global property from a LAMMPS instance.

valid values for lmp_type:resulting return type:
LAMMPS_INTUnsafeArray{Int32, 1}
LAMMPS_DOUBLEUnsafeArray{Float64, 1}
LAMMPS_INT64UnsafeArray{Int64, 1}
LAMMPS_STRINGString

Scalar values get returned as a vector with a single element. This way it's possible to modify the internal state of the LAMMPS instance even if the data is scalar.

Info

Closing the LAMMPS instance or issuing a clear command after calling this method will result in the returned data becoming invalid. To prevent this, copy the returned data.

Warning

Modifying the data through extract_global may lead to inconsistent internal data and thus may cause failures or crashes or bogus simulations. In general it is thus usually better to use a LAMMPS input command that sets or changes these parameters. Those will take care of all side effects and necessary updates of settings derived from such settings.

A full list of global variables can be found in the lammps documentation.

source
LAMMPS.extract_settingMethod
extract_setting(lmp::LMP, name::String)::Int32

Query LAMMPS about global settings.

A full list of settings can be found in the lammps documentation.

Examples

    LMP(["-screen", "none"]) do lmp
        command(lmp, """
            region cell block 0 3 0 3 0 3
            create_box 1 cell
            lattice sc 1
            create_atoms 1 region cell
        """)

        extract_setting(lmp, "dimension") |> println # 3
        extract_setting(lmp, "nlocal") |> println # 27
    end
source
LAMMPS.extract_variableFunction
extract_variable(lmp::LMP, name::String, lmp_variable::LMP_VARIABLE, group::Union{String, Nothing}=nothing)

Extracts the data from a LAMMPS variable. When the variable is either an equal-style compatible variable, a vector-style variable, or an atom-style variable, the variable is evaluated and the corresponding value(s) returned. Variables of style internal are compatible with equal-style variables, if they return a numeric value. For other variable styles, their string value is returned.

valid values for lmp_variable:return type
VAR_ATOMVector{Float64}(copy)
VAR_EQUALFloat64
VAR_STRINGString
VAR_VECTORUnsafeArray{Float64, 1}

the kwarg group determines for which atoms the variable will be extracted. It's only aplicable for VAR_ATOM and will cause an error if used for other variable types. The entires for all atoms not in the group will be zeroed out. By default, all atoms will be extracted.

source
LAMMPS.fix_neighborlistMethod
fix_neighborlist(lmp::LMP, id::String; request=0)

Retrieve neighbor list requested by a fix.

The neighbor list request from a fix is identified by the fix ID and the request ID. The request ID is typically 0, but will be > 0 in case a fix has multiple neighbor list requests.

Each neighbor list contains vectors of local indices of neighboring atoms. These can be used to index into Arrays returned from extract_atom.

source
LAMMPS.gatherFunction
gather(lmp::LMP, name::String, T::Union{Type{Int32}, Type{Float64}}, ids::Union{Nothing, Array{Int32}}=nothing)

Gather the named per-atom, per-atom fix, per-atom compute, or fix property/atom-based entities from all processes. By default (when ids=nothing), this method collects data from all atoms in consecutive order according to their IDs. The optional parameter ids determines for which subset of atoms the requested data will be gathered. The returned data will then be ordered according to ids

Compute entities have the prefix c_, fix entities use the prefix f_, and per-atom entites have no prefix.

The returned Array is decoupled from the internal state of the LAMMPS instance.

ids

The optional parameter ids only works, if there is a map defined. For example by doing: command(lmp, "atom_modify map yes") However, LAMMPS only issues a warning if that's the case, which unfortuately cannot be detected through the underlying API. Starting form LAMMPS version 17 Apr 2024 this should no longer be an issue, as LAMMPS then throws an error instead of a warning.

source
LAMMPS.gather_anglesMethod
gather_angles(lmp::LMP)

Gather the list of all angles into a 4 x nangles Matrix:

row1 -> angle type
row2 -> atom 1
row3 -> atom 2
row4 -> atom 3
source
LAMMPS.gather_bondsMethod
gather_bonds(lmp::LMP)

Gather the list of all bonds into a 3 x nbonds Matrix:

row1 -> bond type
row2 -> atom 1
row3 -> atom 2
source
LAMMPS.gather_dihedralsMethod
gather_dihedrals(lmp::LMP)

Gather the list of all dihedrals into a 5 x ndihedrals Matrix:

row1 -> dihedral type
row2 -> atom 1
row3 -> atom 2
row4 -> atom 3
row5 -> atom 4
source
LAMMPS.gather_impropersMethod
gather_impropers(lmp::LMP)

Gather the list of all impropers into a 5 x nimpropers Matrix:

row1 -> improper type
row2 -> atom 1
row3 -> atom 2
row4 -> atom 3
row5 -> atom 4
source
LAMMPS.get_category_idsFunction
get_category_ids(lmp::LMP, category::String, buffer_size::Integer=50)

Look up the names of entities within a certain category.

Valid categories are: compute, dump, fix, group, molecule, region, and variable. names longer than buffer_size will be truncated to fit inside the buffer.

source
LAMMPS.get_mpi_commMethod
get_mpi_comm(lmp::LMP)::Union{Nothing, MPI.Comm}

Return the MPI communicator used by the lammps instance or nothing if the lammps build doesn't support MPI.

source
LAMMPS.get_natomsMethod
get_natoms(lmp::LMP)::Int64

Get the total number of atoms in the LAMMPS instance.

Will be precise up to 53-bit signed integer due to the underlying lammps_get_natoms returning a Float64.

source
LAMMPS.locateMethod
locate()

Locate the LAMMPS library currently being used, by LAMMPS.jl

source
LAMMPS.pair_neighborlistMethod
pair_neighborlist(lmp::LMP, style::String; exact=false, nsub=0, request=0)

This function determines which of the available neighbor lists for pair styles matches the given conditions. It first matches the style name. If exact is true the name must match exactly, if exact is false, a regular expression or sub-string match is done. If the pair style is hybrid or hybrid/overlay the style is matched against the sub styles instead. If a the same pair style is used multiple times as a sub-style, the nsub argument must be > 0 and represents the nth instance of the sub-style (same as for the pair_coeff command, for example). In that case nsub=0 will not produce a match and this function will Error.

The final condition to be checked is the request ID (reqid). This will normally be 0, but some pair styles request multiple neighbor lists and set the request ID to a value > 0.

Each neighbor list contains vectors of local indices of neighboring atoms. These can be used to index into Arrays returned from extract_atom.

Examples

lmp = LMP()

command(lmp, """
    region cell block 0 3 0 3 0 3
    create_box 1 cell
    lattice sc 1
    create_atoms 1 region cell
    mass 1 1

    pair_style zero 1.0
    pair_coeff * *

    run 1
""")

x = extract_atom(lmp, "x", LAMMPS_DOUBLE_2D; with_ghosts=true)

for (iatom, neighs) in pair_neighborlist(lmp, "zero")
    for jatom in neighs
        ix = @view x[:, iatom]
        jx = @view x[:, jatom]

        println(ix => jx)
   end
end
source
LAMMPS.reset_boxFunction
reset_box(lmp::LMP, boxlo, boxhi, xy::Real = 0, yz::Real = 0, xz::Real = 0)

Reset simulation box parameters.

source
LAMMPS.scatter!Method
scatter!(lmp::LMP, name::String, data::VecOrMat{T}, ids::Union{Nothing, Array{Int32}}=nothing) where T<:Union{Int32, Float64}

Scatter the named per-atom, per-atom fix, per-atom compute, or fix property/atom-based entity in data to all processes. By default (when ids=nothing), this method scatters data to all atoms in consecutive order according to their IDs. The optional parameter ids determines to which subset of atoms the data will be scattered.

Compute entities have the prefix c_, fix entities use the prefix f_, and per-atom entites have no prefix.

ids

The optional parameter ids only works, if there is a map defined. For example by doing: command(lmp, "atom_modify map yes") However, LAMMPS only issues a warning if that's the case, which unfortuately cannot be detected through the underlying API. Starting form LAMMPS version 17 Apr 2024 this should no longer be an issue, as LAMMPS then throws an error instead of a warning.

source
LAMMPS.set_energy!Method
set_energy!(fix::FixExternal, energy::AbstractVector{Float64}; set_global=false)

Sets the contribution to the local per-atom and total global energy resulting from the FixExternal. This function should only be used within the callback of the FixExternal

source
LAMMPS.set_library!Method
set_library!(path)

Change the library path used by LAMMPS.jl for liblammps.so to path.

Note

You will need to restart Julia to use the new library.

Warning

Due to a bug in Julia (until 1.6.5 and 1.7.1), setting preferences in transitive dependencies is broken (https://github.com/JuliaPackaging/Preferences.jl/issues/24). To fix this either update your version of Julia, or add LAMMPS_jll as a direct dependency to your project.

source
LAMMPS.set_virial!Method
set_virial!(fix::FixExternal, virial_peratom::AbstractVector{SVector{6, Float64}})
set_virial!(fix::FixExternal, virial_peratom::AbstractMatrix{Float64})

Sets the contribution to the local per-atom and total global virial stress tensors resulting from the FixExternal. This function should only be used within the callback of the FixExternal

The virial tensor is represented as a 6-component vector with the order: [xx, yy, zz, xy, xz, yz].

source