API

LAMMPS.LMPType
LMP(f::Function, args=String[], comm=nothing)

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::Union{Nothing, MPI.Comm}=nothing)

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.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.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; copy=false, with_ghosts=false)

Extract per-atom data from the lammps instance.

valid values for lmp_type:resulting return type:
LAMMPS_INTVector{Int32}
LAMMPS_INT_2DMatrix{Int32}
LAMMPS_DOUBLEVector{Float64}
LAMMPS_DOUBLE_2DMatrix{Float64}
LAMMPS_INT64Vector{Int64}
LAMMPS_INT64_2DMatrix{Int64}
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, set copy=true.

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

Arguments

  • copy: determines whether lammps internal memory is used or if a copy is made.
  • with_ghosts: Determines wheter entries for ghost atoms are included. This is ignored for "mass".
source
LAMMPS.extract_computeMethod
extract_compute(lmp::LMP, name::String, style::_LMP_STYLE_CONST, lmp_type::_LMP_TYPE; copy::Bool=true)

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_SCALARVector{Float64}
TYPE_VECTORVector{Float64}
TYPE_ARRAYMatrix{Float64}
SIZE_VECTORVector{Int32}
SIZE_COLSVector{Int32}
SIZE_ROWSVector{Int32}

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

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, set copy=true.

Examples

LMP(["-screen", "none"]) do lmp
    extract_compute(lmp, "thermo_temp", LMP_STYLE_GLOBAL, TYPE_VECTOR, copy=true)[2] = 2
    extract_compute(lmp, "thermo_temp", LMP_STYLE_GLOBAL, TYPE_VECTOR, copy=false)[3] = 3

    extract_compute(lmp, "thermo_temp", LMP_STYLE_GLOBAL, TYPE_SCALAR) |> println # [0.0]
    extract_compute(lmp, "thermo_temp", LMP_STYLE_GLOBAL, TYPE_VECTOR) |> println # [0.0, 0.0, 3.0, 0.0, 0.0, 0.0]
end
source
LAMMPS.extract_globalMethod
extract_global(lmp::LMP, name::String, lmp_type::_LMP_DATATYPE; copy::Bool=false)

Extract a global property from a LAMMPS instance.

valid values for lmp_type:resulting return type:
LAMMPS_INTVector{Int32}
LAMMPS_INT_2DMatrix{Int32}
LAMMPS_DOUBLEVector{Float64}
LAMMPS_DOUBLE_2DMatrix{Float64}
LAMMPS_INT64Vector{Int64}
LAMMPS_INT64_2DMatrix{Int64}
LAMMPS_STRINGString (always a copy)

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, set copy=true.

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; copy::Bool=false)

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}
VAR_EQUALFloat64
VAR_STRINGString
VAR_VECTORVector{Float64}

the kwarg copy determies wheter a copy of the underlying data is made. copy is only aplicable for VAR_VECTOR and VAR_ATOM. For all other variable types, a copy will be made regardless. The underlying LAMMPS API call for VAR_ATOM internally allways creates a copy of the data. As the memory for this gets allocated by LAMMPS instead of julia, it needs to be dereferenced using LAMMPS.API.lammps_free instead of through the garbage collector. If copy=false this gets acieved by registering LAMMPS.API.lammps_free as a finalizer for the returned data. Alternatively, setting copy=true will instead create a new copy of the data. The lammps allocated block of memory will then be freed immediately.

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.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_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.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_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