ACEpotentials.jl API

Exported

ACEpotentials.acefit!Method

function acefit!(model, data; kwargs...) : provides a simplified interface to fitting the parameters of a model specified via ACE1Model. The data should be provided as a collection (AbstractVector) of JuLIP.Atoms structures.

Keyword arguments:

  • energy_key, force_key, virial_key specify

the label of the data to which the parameters will be fitted.

  • weights specifies the regression weights, default is 30 for energy, 1 for forces and virials
  • solver specifies the lsq solver, default is BLR (BayesianLinearRegression)
  • smoothness specifies the smoothness prior, i.e. how strongly damped parameters corresponding to high polynomial degrees are; is 2.
  • prior specifies a covariance of the prior, if nothing then a smoothness prior is used, using the smoothness parameter
  • repulsion_restraint specifies whether to add artificial data to the training set that effectively introduces a restraints encouraging repulsion in the limit rij -> 0.
  • restraint_weight specifies the weight of the repulsion restraint.
  • export_lammps : path to a file to which the fitted potential will be exported in a LAMMPS compatible format (yace)
  • export_json : path to a file to which the fitted potential will be exported in a JSON format, which can be read from Julia or Python
source
ACEpotentials.load_potentialMethod
function load_potential(fname::AbstractString; new_format=false, verbose=true)

Load ACE potential from given file fname.

Kwargs

  • new_format=false - If true returns potential as ACEmd.ACEpotential format, else use old JuLIP format
  • verbose=true - Display version info on load
source
ACEpotentials.save_potentialMethod
save_potential( fname, potential::ACE1x.ACE1Model; save_version_numbers=true, meta=nothing)

Save ACE potentials. Prefix is either .json, .yml or .yace, which also determines file format.

Kwargs

  • saveversionnumbers=true : If true save version information or relevant packages
  • meta=nothing : Seve some metadata with the potential (needs to be Dict{String, Any})
source
ACEpotentials.site_descriptorsFunction
site_descriptors(basis, atoms::AbstractAtoms[, domain])

Compute site descriptors for all atoms in atoms, returning them as a vector of vectors. If the optional argument domain is passed as a list of integers (atom indices), then only the site descriptors for those atoms are computed and returned.

source

Not exported

ACEpotentials.PolyTransform_paramsMethod

PolyTransform_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct ACE1.Transforms.PolyTransform. All parameters are passed as keyword argument. Also see?PolyTransform`

Implements the distance transform

\[ x(r) = \Big(\frac{1 + r_0}{1 + r}\Big)^p\]

Parameters

  • p = 2
  • r0 = 2.5
source
ACEpotentials.ace_basis_paramsMethod

ace_basis_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct an ACE basis (RPIBasis). All parameters are passed as keyword argument. If no default is given then the argument is required.

Parameters

  • species : single species or list of species (mandatory)
  • N : correlation order, positive integer (mandatory)
  • maxdeg : maximum degree, positive real number (note the precise notion of degree is specified by further parameters) (mandatory)
  • r0 = 2.5 : rough estimate for nearest neighbour distance
  • radial = radial_basis_params(; r0 = r0) : one-particle basis parameters; cf ?radial_basis_params for details
  • transform = transform_params(; r0 = r0) : distance transform parameters; cf ?transform_params() for details
  • degree = degree_params() : class of sparse polynomial degree to select the basis; see ?degree_params for details
source
ACEpotentials.at_dimerMethod

function at_dimer(r, z1, z0) : generates a dimer with separation r and atomic numbers z1 and z0. (can also use symbols or strings)

source
ACEpotentials.at_trimerMethod

function at_trimer(r1, r2, θ, z0, z1, z2) : generates a trimer with separations r1 and r2, angle θ and atomic numbers z0, z1 and z2 (can also use symbols or strings), where z0 is the species of the central atom, z1 at distance r1 and z2 at distance r2.

source
ACEpotentials.basis_paramsMethod

basis_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct one of the basis. All parameters are passed as keyword argument and the kind of parameters required depend on "type".

ACE (RPI) basis

Returns a dictionary containing the complete set of parameters required to construct an ACE basis (RPIBasis). All parameters are passed as keyword argument. If no default is given then the argument is required.

Parameters

  • type = "ace"
  • species : single species or list of species (mandatory)
  • N : correlation order, positive integer (mandatory)
  • maxdeg : maximum degree, positive real number (note the precise

notion of degree is specified by further parameters) (mandatory)

  • r0 = 2.5 : rough estimate for nearest neighbour distance
  • radial = radial_basis_params(; r0 = r0) : one-particle basis

parameters; cf ?basis_params of type "radial" for details

  • transform = transform_params(; r0 = r0) : distance transform

parameters; cf ?transform_params() for details

  • degree = degree_params() : class of sparse polynomial degree

to select the basis; see ?degree_params for details

Pair basis

Returns a dictionary containing the complete set of parameters required to construct an pair basis (PolyPairBasis). All parameters are passed as keyword argument.

Parameters

  • type = "pair"
  • species : single species or list of species (mandatory)
  • maxdeg : maximum degree, positive real number (note the precise

notion of degree is specified by further parameters) (mandatory)

  • r0 = 2.5 : rough estimate for nearest neighbour distance
  • rcut = 5.0: outer cutoff, Å
  • rin = 0.0: inner cutoff, Å
  • pcut = 2: outer cutoff parameter; * pcut=2: function and first derivative go to zero at the outer cutoff * pcut=1: function forced to go through zero at the outer cutoff * pcut=0: no constraint at the outer cutoff
  • pin = 0: inner cutoff parameter * pin=2: function and first derivative go to zero at the inner cutoff * pin=1: function forced to go through zero at the inner cutoff * pin=0: no constraint at the inner cutoff
  • transform = transform_params(; r0 = r0) : distance transform

parameters; cf ?transform_params() for details

Radial basis of ACE

Returns a dictionary containing the complete set of parameters required to construct radial basis for ACE. All parameters are passed as keyword argument.

Parameters

  • type = "radial"
  • r0 = 2.5 : rough estimate for nearest neighbour distance
  • rcut = 5.0: outer cutoff, Å
  • rin = 0.5 * r0: inner cutoff, Å
  • pcut = 2: outer cutoff parameter; * pcut=2: function and first derivative go to zero at the outer cutoff * pcut=1: function forced to go through zero at the outer cutoff * pcut=0: no constraint at the outer cutoff
  • pin = 2: inner cutoff parameter * pin=2: function and first derivative go to zero at the inner cutoff * pin=1: function forced to go through zero at the inner cutoff * pin=0: no constraint at the inner cutoff
source
ACEpotentials.blr_paramsMethod

blr_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a Bayesian linear regression solver. All parameters are passed as keyword argument.

Parameters

  • verbose = false
source
ACEpotentials.data_paramsMethod

data_params(; kwargs...)` : returns a dictionary containing the complete set of parameters required to read data from an .xyz file. j All parameters are passed as keyword arguments.

Parameters

  • fname : a "*.xyz" file with atomistic data (mandatory).
  • energy_key = "energy" : key identifying energy for fitting.
  • force_key = "forces" : key identifying forces for fitting.
  • virial_key = "virial" : key identifying virial tensor for fitting.
  • weight_key = "config_type` : key identifying label for setting the correct weight from weights dictionary.
source
ACEpotentials.db_paramsMethod

db_params(; kwargs...)` : returns a dictionary containing all of the parameters needed for making a LsqDB. All parameters are passed as keyword argumts.

Parameters

  • data : data parameters, see ?data_params for details (mandatory)

  • basis : dictionary containing dictionaries that specify the basis used in fitting. For example

    julia basis = Dict( "pair_short" => Dict( "type" => "pair", ...), "pair_long" => Dict("type" => "pair", ...), "manybody" => Dict("type" => "ace", ...), "nospecies" => Dict("type" => "ace", species = ["X",], ...)

keys of basis are ignored, so that multiple basis with different specifications (e.g. smaller and larger cutoffs) can be combined. See ?basis_params for more detail.

  • LSQ_DB_fname_stem = "" : stem to save LsqDB to. Doesn't get saved if set to an empty string (""). If LSQ_DB_fname_stem * "_kron.h5" file is not present it gets renamed, a new LsqDB is constructed and saved.
source
ACEpotentials.decohesion_curveMethod

Generate a decohesion curve for testing the smoothness of a potential. Arguments:

  • at0 : unit cell
  • pot : potential implementing energy

Keyword Arguments:

  • dim = 1 : dimension into which to expand
  • mult = 10 : multiplicative factor for expanding the cell in dim direction
  • aa = :auto : array of stretch values of the lattice parameter to use
  • npoints = 100 : number of points to use in the stretch array (for auto aa)
source
ACEpotentials.degree_paramsMethod

degree_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a specification for polynomial degree. All parameters are passed as keyword argument and the kind of parameters required depend on "type".

SparsePSHDegree

Returns a dictionary containing the complete set of parameters required to construct ACE1.RPI.SparsePSHDegree. See ?SparsePSHDegree.

Parameters

  • type = "sparse"
  • wL = 1.5
  • csp = 1.0
  • chc = 0.0
  • chc = 0.0
  • ahc = 0.0
  • bhc = 0.0
  • p = 1.0

##SparsePSHDegreeM Returns a dictionary containing the complete set of parameters required to construct ACE1.RPI.SparsePSHDegree. Also see ?SparsePSHDegreeM.

NB maxdeg of ACE basis (RPIBasis) has to be set to 1.0.

Parameters

  • Dd : Dictionary specifying max degrees (mandatory)
  • Dn = Dict("default" => 1.0) : Dictionary specifying weights for degree

of radial basis functions (n)

  • Dl = Dict("default" => 1.5) : Dictionary specifying weights for degree

of angular basis functions (l)

Each dictionary should have a "default" entry. In addition, different degrees or weights can be specified for each correlation order and/or correlation order-species combination. For example

"Dd" => Dict(
      "default" => 10,
      3 => 9,
      (4, "C") => 8,
      (4, "H") => 0)

in combination with N=4 and maxdeg=1.0, will set maximum polyonmial degree on N=1 and N=2 functions to 10, to 9 for N=3 functions and will only allow N=4 basis functions on carbon atoms, up to polynomial degree 8.

source
ACEpotentials.dimer_energyMethod

function dimer_energy(pot, r, z1, z0) : energy of a dimer with separation r and atomic numbers z1 and z0 using the potential pot; subtracting the 1-body contributions.

source
ACEpotentials.dimersMethod

dimers(potential, elements; kwargs...) : Generate a dictionary of dimer curves for a given potential.

  • potential : potential to use to evaluate energy
  • elements : list of chemical species, symbols for which the dimers are to be computed

The function returns a dictionary Ddim such that D[(s1, s2)] contains pairs or arrays (rr, E) which can be plotted plot(rr, E).

source
ACEpotentials.fill_defaultsMethod

fill_defaults(params::Dict; param_key = "fit_params") -> params

Recursively updates any missing entries with default parameters. Accepted param_key values and corresponding functions:

    "fit_params" => ACEpotentials.fit_params,
    "data" => ACEpotentials.data_params,
    "solver" => ACEpotentials.solver_params,
    "basis" => ACEpotentials.basis_params,
    "ace" => ACEpotentials.ace_basis_params,  
    "pair" => ACEpotentials.pair_basis_params,
    "radial" => ACEpotentials.radial_basis_params,
    "transform" => ACEpotentials.transform_params, 
    "degree" => ACEpotentials.degree_params,
    "P" => ACEpotentials.regularizer_params
source
ACEpotentials.fit_aceMethod

fit_ace(params::Dict) -> IP, lsqinfo

Function to set up and fit the least-squares problem of "atoms' positions" -> "energy, forces, virials". Takes in a dictionary with all the parameters. See ?fit_params for details.

source
ACEpotentials.fit_paramsMethod

fit_params(; kwargs...)

Returns a dictionary containing all of the parameters needed to make an ACE potential. All parameters are passed as keyword argumts.

Parameters

  • data : data parameters, see ?data_params for details (mandatory)
  • basis : dictionary containing dictionaries that specify the basis used in fitting. For example
basis = Dict(
    "pair_short" => Dict( "type" => "pair", ...), 
    "pair_long" => Dict("type" => "pair", ...), 
    "manybody" => Dict("type" => "ace", ...), 
    "nospecies" => Dict("type" => "ace", species = ["X",], ...)

keys of basis are ignored, so that multiple basis with different specifications (e.g. smaller and larger cutoffs) can be combined. See ?basis_params for more detail.

  • solver : dictionary containing parameters that specify the solver for least squares problem (mandatory). See ?solver_params.
  • e0 : Dict{String, Float} containing reference values for isolated atoms' energies (mandatory).
  • weights : dictionary of Dict("config_type" => Dict("E" => Float, "F => Float))entries specifying fitting weights. "default" is set to1.0` for all of "E", "F", and "V" weights.
  • P : regularizer parameters (optional), see ?regularizer_params.
  • ACE_fname = "ACE_fit.json" : filename to save ACE to. Potential & info do not get saved if ACE_fname isnothing() or is set to "". Files already parseentry are renamed and not overwritten.
  • LSQ_DB_fname_stem = "" : stem to save LsqDB to. Doesn't get saved if set to an empty string (""). If the file is already present, but fit_from_LSQ_DB is set to false, the old database is renamed, a new one constructed and saved under the given name.
  • fit_from_LSQ_DB = false: whether to fit from a least squares database specified with LSQ_DB_fname_stem. If LSQ_DB_fname_stem * "_kron.h5" file is not present, LsqDB is constructed from scratch and saved.
source
ACEpotentials.get_adfMethod

function get_adf(data::AbstractVector{<: Atoms}, r_cut; kwargs...) :

Angular distribution, i.e. list of angles in [0, π] between all pairs of bonds of length at most r_cut. Keyword arguments:

  • skip = 3 : only consider every skipth atom in the dataset.
  • maxsamples = 100_000 : maximum number of samples to return.
source
ACEpotentials.get_rdfMethod

function get_rdf(data::AbstractVector{<: Atoms}, r_cut; kwargs...) :

Produce a list of r values that occur in the dataset, restricted to the cutoff radius r_cut. Keyword arguments:

  • rescale = true : resample the data to account for volume scaling, i.e. a distance r will be kept with probability min(1, (r0/r)^2).
  • r0 = :min : parameter for resampling. If :min then the minimum r occuring in the dataset is taken.
  • maxsamples = 100_000 : maximum number of samples to return.
source
ACEpotentials.laplacian_regularizer_paramsMethod

laplacian_regularizer_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a laplacian regularizer. All parameters are passed as keyword argument.

Parameters

  • rlap_scal = 3.0
source
ACEpotentials.lsqr_paramsMethod

lsqr_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a lsqr solver. All parameters are passed as keyword argument.

Parameters

  • damp = 5e-3
  • atol = 1e-6
  • colim = 1e8
  • maxiter = 1e5
  • verbose = false
source
ACEpotentials.make_ace_dbMethod

make_ace_db(params::Dict) -> LsqDB

Makes a LsqDB from given parameters' dictionary. For params see ?db_params; parameters from fit_params also work, except unnecessary entries will be ignored. Returns IPFitting.LsqDB

source
ACEpotentials.multitransform_paramsMethod

multitransform_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct ACE.Transform.multitransform. All parameters are passed as keyword argument.

Parameters

  • transforms : dictionary specifying transforms for each species pair. Can be

given per-pair (i.e. only for ("element1", "element2") and not for ("element2", "element1")) or can be different for ("element1", "element2") and ("element2", "element1"). For example

transforms = Dict(
      ("C", "C") => Dict("type"=> "polynomial"),
      ("C", "H") => Dict("type"=> "polynomial"),
      ("H", "H") => Dict("type" => "polynomial"))
  • rin, rcut: values for inner and outer cutoffs, alternative to cutoffs
  • cutoffs : dictionary specifying inner and outer cutoffs for each element pair

(either symmetrically or non-symmetrically). Alternative to rin & rcut. For example

cutoffs => Dict(
      ("C", "C") => (1.1, 4.5),
      ("C", "H") => (0.9, 4.5),
      ("H", "H") => (1.23, 4.5)),
source
ACEpotentials.pair_basis_paramsMethod

pair_basis_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct an pair basis (PolyPairBasis). All parameters are passed as keyword argument.

Parameters

  • species : single species or list of species (mandatory)
  • maxdeg : maximum degree, positive real number (note the precise notion of degree

is specified by further parameters) (mandatory)

  • r0 = 2.5 : rough estimate for nearest neighbour distance
  • rcut = 5.0: outer cutoff, Å
  • rin = 0.0: inner cutoff, Å
  • pcut = 2: outer cutoff parameter; * pcut=2: function and first derivative go to zero at the outer cutoff * pcut=1: function forced to go through zero at the outer cutoff * pcut=0: no constraint at the outer cutoff
  • pin = 0: inner cutoff parameter * pin=2: function and first derivative go to zero at the inner cutoff * pin=1: function forced to go through zero at the inner cutoff * pin=0: no constraint at the inner cutoff
  • transform = transform_params(; r0 = r0) : distance transform parameters;

cf ?transform_params() for details

source
ACEpotentials.parse_ace_basis_keysMethod

parse_ace_basis_keys(ace_basis::Dict) -> ace_basis

("C", "C")-type tuples are saved to and read back in from JSON as "("C", "C")" .json. It's slightly easier to save these to JSON or YAM as "(C, C)". This function converts "(C, C)"-type strings back to parameter-friendly ("C", "C").

source
ACEpotentials.qr_paramsMethod

qr_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a qr solver. All parameters are passed as keyword argument.

source
ACEpotentials.radial_basis_paramsMethod

radial_basis_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct radial basis for ACE. All parameters are passed as keyword argument.

Parameters

  • r0 = 2.5 : rough estimate for nearest neighbour distance
  • rcut = 5.0: outer cutoff, Å
  • rin = 0.5 * r0: inner cutoff, Å
  • pcut = 2: outer cutoff parameter; * pcut=2: function and first derivative go to zero at the outer cutoff * pcut=1: function forced to go through zero at the outer cutoff * pcut=0: no constraint at the outer cutoff
  • pin = 2: inner cutoff parameter * pin=2: function and first derivative go to zero at the inner cutoff * pin=1: function forced to go through zero at the inner cutoff * pin=0: no constraint at the inner cutoff
source
ACEpotentials.regularizer_paramsMethod

regularizer_params(; type = "laplacian", kwargs...) : returns a dictionary containing the complete set of parameters required to construct one of the solvers. All parameters are passed as keyword argument and the kind of parameters required depend on "type".

LSQR Parameters (default)

  • type = "laplacian"
  • rlap_scal = 3.0
source
ACEpotentials.rrqr_paramsMethod

rrqr_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a rrqr solver. All parameters are passed as keyword argument.

Parameters

  • rtol = 1e-5
source
ACEpotentials.save_fitMethod

save_fit(fname, IP, lsqinfo)

Saves Dict("IP" => IP, "info" => lsqinfo) to fname. If fname is already present, it is renamed and dictionary saved to fname.

source
ACEpotentials.sklearn_ard_paramsMethod

sklearn_ard_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a scikit-learn Automatic Relevance Detemrination solver. All parameters are passed as keyword argument.

Parameters

  • n_iter = 300
  • tol = 1e-3
  • threshold_lambda = 1e4
source
ACEpotentials.sklearn_brr_paramsMethod

sklearn_brr_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct a scikit-learn Bayesian ridge regression solver. All parameters are passed as keyword argument.

Parameters

  • n_iter = 300
  • tol = 1e-3
source
ACEpotentials.solver_paramsMethod

solver_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct one of the solvers. All parameters are passed as keyword argument and the kind of parameters required depend on "type".

QR Parameters

  • 'type = "qr"`

LSQR Parameters

  • type = "lsqr"
  • damp = 5e-3
  • atol = 1e-6
  • colim = 1e8
  • maxiter = 1e5
  • verbose = false

RRQR Parameters

  • type = "rrqr"
  • rtol = 1e-5

SKLEARN_BRR

  • type = "sklearn_brr"
  • n_iter = 300
  • tol = 1e-3

SKLEARN_ARD

  • type = "sklearn_ard"
  • n_iter = 300
  • tol = 1e-3
  • threshold_lambda = 1e4

BLR

  • type = "blr"
  • verbose = false
source
ACEpotentials.sparse_degree_M_paramsMethod

sparse_degree_M_params(;kwargs...): Returns a dictionary containing the complete set of parameters required to construct ACE1.RPI.SparsePSHDegree. Also see ?SparsePSHDegreeM.

NB maxdeg of ACE basis (RPIBasis) has to be set to 1.0.

Parameters

  • Dd : Dictionary specifying max degrees (mandatory)
  • Dn = Dict("default" => 1.0) : Dictionary specifying weights for degree

of radial basis functions (n)

  • Dl = Dict("default" => 1.5) : Dictionary specifying weights for degree

of angular basis functions (l)

Each dictionary should have a "default" entry. In addition, different degrees or weights can be specified for each correlation order and/or correlation order-species combination. For example

"Dd" => Dict(
      "default" => 10,
      3 => 9,
      (4, "C") => 8,
      (4, "H") => 0)

in combination with N=4 and maxdeg=1.0, will set maximum polyonmial degree on N=1 and N=2 functions to 10, to 9 for N=3 functions and will only allow N=4 basis functions on carbon atoms, up to polynomial degree 8.

source
ACEpotentials.sparse_degree_paramsMethod

sparse_degree_params(; kwargs...): returns a dictionary containing the complete set of parameters required to construct ACE1.RPI.SparsePSHDegree. See ?SparsePSHDegree.

Parameters

  • wL = 1.5
  • csp = 1.0
  • chc = 0.0
  • chc = 0.0
  • ahc = 0.0
  • bhc = 0.0
  • p = 1.0

NB p = 1 is current ignored, but we put it in so we can experiment later with p = 2, p = inf.

source
ACEpotentials.transform_paramsMethod

transform_params(; kwargs...) : returns a dictionary containing the complete set of parameters required to construct one of the transforms. All parameters are passed as keyword argument and the kind of parameters required depend on "type".

Polynomial transform

Returns a dictionary containing the complete set of parameters required to construct ACE1.Transforms.PolyTransform. All parameters are passed as keyword argument. Also see?PolyTransform`

Implements the distance transform

\[ x(r) = \Big(\frac{1 + r_0}{1 + r}\Big)^p\]

Parameters

  • type = "polynomial"
  • p = 2
  • r0 = 2.5

Multitransform

Returns a dictionary containing the complete set of parameters required to construct ACE.Transform.multitransform. All parameters are passed as keyword argument.

Parameters

  • transforms : dictionary specifying transforms for each species pair. Can be

given per-pair (i.e. only for ("element1", "element2") and not for ("element2", "element1")) or can be different for ("element1", "element2") and ("element2", "element1"). For example

transforms = Dict(
      ("C", "C") => Dict("type"=> "polynomial"),
      ("C", "H") => Dict("type"=> "polynomial"),
      ("H", "H") => Dict("type" => "polynomial"))
  • rin, rcut: values for inner and outer cutoffs, alternative to cutoffs
  • cutoffs : dictionary specifying inner and outer cutoffs for each element pair

(either symmetrically or non-symmetrically). Alternative to rin & rcut. For example

cutoffs => Dict(
      ("C", "C") => (1.1, 4.5),
      ("C", "H") => (0.9, 4.5),
      ("H", "H") => (1.23, 4.5)),

identity

IdTransform_params(;) : returns Dict("type" => "identity"), needed to construct ACE1.Transforms.IdTransform.

source
ACEpotentials.trimer_energyMethod

function trimer_energy(IP, r1, r2, θ, z0, z1, z2) : computes the energy of a trimer, subtracting the 2-body and 1-body contributions.

source
ACEpotentials.trimersMethod

trimers(potential, elements, r1, r2; kwargs...) : Generate a dictionary of trimer curves for a given potential.

  • potential : potential to use to evaluate energy
  • elements : list of chemical species, symbols for which the trimers are to be computed
  • r1, r2 : distance between the central atom and the first, second neighbour

The function returns a dictionary Dtri such that D[(s1, s2, s3)] contains pairs or arrays (θ, E) which can be plotted plot(θ, E).

source