API reference
GeometryOptimization.Autoselect
GeometryOptimization.DofManager
GeometryOptimization.GeoOptDefaultCallback
GeometryOptimization.GeoOptProblem
GeometryOptimization.OptimCG
GeometryOptimization.OptimLBFGS
GeometryOptimization.OptimSD
GeometryOptimization.analyze_mask
GeometryOptimization.clamp_atoms
GeometryOptimization.convert_to_updatable
GeometryOptimization.minimize_energy!
GeometryOptimization.Autoselect
— TypeUse a heuristic to automatically select the minimisation algorithm (Currently OptimCG
, but this may change silently)
GeometryOptimization.DofManager
— TypeDofManager
:
Constructor:
DofManager(sys; variablecell = false, r0 =..., free=..., clamp=..., mask=...)
variablecell
determines whether the cell is fixed or allowed to change during optimizationr0
is a reference length-scale, default set to one (in the unit of sys), this is used to non-dimensionalize the degrees of freedom.
In addition set at most one of the kwargs:
- no kwarg: all atoms are free
free
: list of free atom indices (not dof indices)clamp
: list of clamped atom indices (not dof indices)mask
: 3 x N Bool array to specify individual coordinates to be clamped
Meaning of dofs
On call to the constructor, DofManager
stores positions and cell X0, C0
, dofs are understood relative to this initial configuration. get_dofs(sys, dm::DofManager)
returns a vector that represents the non-dimensional displacement and a deformation matrix (U, F)
. The new configuration extracted from a dof vector is understood as
- The new cell:
C = F * C0
- The new positions:
𝐫[i] = F * (X0[i] + U[i] * r0)
One aspect of this definition is that clamped atom positions still change via the deformation F
. This is natural in the context of optimizing the cell shape.
GeometryOptimization.GeoOptDefaultCallback
— TypeCallback producing a convergence table summarising the geometry optimisation convergence. If always_show_header=true
the header is shown in each iteration. This is helpful if the calculator produces output as well.
GeometryOptimization.GeoOptProblem
— TypeGeoOptProblem(system, calculator, dofmgr, geoopt_state; kwargs...)
Turn system
, calculator
and geoopt_state::GeometryOptimizationState
into an OptimizationProblem
for solve_problem
. Note that the system
is not updated automatically and that internally atomic units are used.
GeometryOptimization.OptimCG
— TypeUse Optim's ConjugateGradient implementation with some good defaults.
GeometryOptimization.OptimLBFGS
— TypeUse Optim's LBFGS implementation with some good defaults.
GeometryOptimization.OptimSD
— TypeUse Optim's GradientDescent (Steepest Descent) implementation with some good defaults.
GeometryOptimization.analyze_mask
— Methodanalyze_mask(sys, free, clamp, mask) -> Any
analyze_mask
: helper function to generate list of dof indices from lists of atom indices indicating free and clamped atoms
GeometryOptimization.clamp_atoms
— Methodclamp_atoms(
system,
clamped_indexes::Union{Nothing, AbstractVector{<:Integer}}
) -> AtomsBase.FlexibleSystem{_A, _B, AtomsBase.PeriodicCell{_A1, _B1}} where {_A, _B, _A1, _B1}
Clamp given atoms in the system. Clamped atoms are fixed and their positions will not be optimized. The atoms to be clamped should be given as a list of indices corresponding to their positions in the system.
This is this a very experimental interface and will likely change in the future.
GeometryOptimization.convert_to_updatable
— Methodconvert_to_updatable(system) -> Any
Convert the input system to a kind of system we can work with, i.e. one where we have the particles
keyword argument for setting new particles and positions.
GeometryOptimization.minimize_energy!
— Functionminimize_energy!(system, calculator; ...)
minimize_energy!(system, calculator, solver; kwargs...)
Minimise the energy of a system
using the specified calculator
. Optimises either only atomic positions (variablecell=false
) or both positions and lattice (variablecell=true
). Returns a named tuple including the optimised system as first entry. Typical arguments passed as solver
are Autoselect()
(the default), OptimLBFGS()
, OptimCG()
, OptimSD()
. These automatically choose some heuristics for setting up the solvers, which we found to work well in practice.
Beyond that any other first-order solver
compatible with Optimization.jl can also be employed. Note, that in principle all such solvers should work, but we only tested a small fraction and you can expect that minor modifications are needed to make some solvers work (PRs appreciated!). In general only first-order or second-order methods work.
See GeometryOptimization.clamp_atoms
for an experimental way how to clamp atoms.
Keyword arguments:
variablecell
: Determines whether the cell is fixed or allowed to change during optimizationmaxiters
: Maximal number of iterationsmaxtime
: Maximal allowed runtime (in seconds)tol_energy
: Tolerance in the energy to stop the minimisation (alltol_*
need to be satisfied)tol_forces
: Tolerance in the force to stop the minimisation (alltol_*
need to be satisfied)tol_virial
: Tolerance in the virial to stop the minimisation (alltol_*
need to be satisfied)maxstep
: Maximal step size (in AU or length units) to be taken in a single optimisation step (not supported for allsolver
s)verbosity
: Printing level. The idea is that0
is silent,1
displays the optimisation progress and≥ 2
starts displaying things from the calculator as well (e.g SCF iterations).callback
: A custom callback, which obtains the pair(optimization_state, geoopt_state)
and is expected to returnfalse
(continue iterating) ortrue
(halt iterations). Note that specifying this overwrites the default printing callback. The calculation thus becomes silent unless aGeoOptDefaultCallback
is included in the callback.kwargs
: All other keyword arguments are passed to the call tosolve
. Note, that if specialkwargs
should be passed to theOptimization.OptimizationProblem
the user needs to setup the problem manually (e.g.OptimizationProblem(system, calculator)
)