Files
lammps/doc/src/min_spin.txt
julient31 f1c3b9d0bf Commit2 JT 072319
- corrected some mistakes in doc files
- modified oso examples to match new line options
2019-07-23 11:24:52 -06:00

97 lines
2.9 KiB
Plaintext

"LAMMPS WWW Page"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c
:link(lws,http://lammps.sandia.gov)
:link(ld,Manual.html)
:link(lc,Commands_all.html)
:line
min_style spin command :h3
min_style spin_oso_cg command :h3
min_style spin_oso_lbfgs command :h3
[Syntax:]
min_style spin
min_style spin_oso_cg
min_style spin_oso_lbfgs :pre
[Examples:]
min_style spin_oso_lbfgs
min_modify discrete_factor 10.0 line spin_none :pre
[Description:]
Apply a minimization algorithm to use when a "minimize"_minimize.html
command is performed.
Style {spin} defines a damped spin dynamics with an adaptive
timestep, according to:
:c,image(Eqs/min_spin_damping.jpg)
with lambda a damping coefficient (similar to a Gilbert
damping).
Lambda can be defined by setting the {alpha_damp} keyword with the
"min_modify"_min_modify.html command.
The minimization procedure solves this equation using an
adaptive timestep. The value of this timestep is defined
by the largest precession frequency that has to be solved in the
system:
:c,image(Eqs/min_spin_timestep.jpg)
with {|omega|_{max}} the norm of the largest precession frequency
in the system (across all processes, and across all replicas if a
spin/neb calculation is performed).
Kappa defines a discretization factor {discrete_factor} for the
definition of this timestep.
{discrete_factor} can be defined with the "min_modify"_min_modify.html
command.
Style {spin_oso_cg} defines an orthogonal spin optimization
(OSO) combined to a conjugate gradient (CG) algorithm.
The "min_modify"_min_modify.html command can be used to
couple the {spin_oso_cg} to a line search procedure, and to modify the
discretization factor {discrete_factor}.
Style {spin_oso_lbfgs} defines an orthogonal spin optimization
(OSO) combined to a limited-memory Broyden-Fletcher-Goldfarb-Shanno
(LBFGS) algorithm.
By default, style {spin_oso_lbfgs} uses a line search procedure.
The "min_modify"_min_modify.html command can be used to
deactivate the line search procedure, and to modify the
discretization factor {discrete_factor}.
For more information about styles {spin_oso_cg} and {spin_oso_lbfgs},
see their implementation reported in "(Ivanov)"_#Ivanov1.
NOTE: All the {spin} styles replace the force tolerance by a torque
tolerance. See "minimize"_minimize.html for more explanation.
NOTE: The {spin_oso_cg} and {spin_oso_lbfgs} styles can be used
for magnetic NEB calculations only if the line search procedure
is deactivated. See "neb/spin"_neb_spin.html for more explanation.
[Restrictions:]
This minimization procedure is only applied to spin degrees of
freedom for a frozen lattice configuration.
[Related commands:]
"min_style"_min_style.html, "minimize"_minimize.html,
"min_modify"_min_modify.html
[Default:]
The option defaults are {alpha_damp} = 1.0 and {discrete_factor} =
10.0.
:line
:link(Ivanov1)
[(Ivanov)] Ivanov, Uzdin, Jonsson. arXiv preprint arXiv:1904.02669, (2019).