2016-05-10 01:33:12 +08:00
<!DOCTYPE html>
<!-- [if IE 8]><html class="no - js lt - ie9" lang="en" > <![endif] -->
<!-- [if gt IE 8]><! --> < html class = "no-js" lang = "en" > <!-- <![endif] -->
< head >
< meta charset = "utf-8" >
< meta name = "viewport" content = "width=device-width, initial-scale=1.0" >
< title > fix ti/rs command — LAMMPS documentation< / title >
< link rel = "stylesheet" href = "_static/css/theme.css" type = "text/css" / >
< link rel = "stylesheet" href = "_static/sphinxcontrib-images/LightBox2/lightbox2/css/lightbox.css" type = "text/css" / >
< link rel = "top" title = "LAMMPS documentation" href = "index.html" / >
< script src = "_static/js/modernizr.min.js" > < / script >
< / head >
< body class = "wy-body-for-nav" role = "document" >
< div class = "wy-grid-for-nav" >
< nav data-toggle = "wy-nav-shift" class = "wy-nav-side" >
< div class = "wy-side-nav-search" >
< a href = "Manual.html" class = "icon icon-home" > LAMMPS
< / a >
< div role = "search" >
< form id = "rtd-search-form" class = "wy-form" action = "search.html" method = "get" >
< input type = "text" name = "q" placeholder = "Search docs" / >
< input type = "hidden" name = "check_keywords" value = "yes" / >
< input type = "hidden" name = "area" value = "default" / >
< / form >
< / div >
< / div >
< div class = "wy-menu wy-menu-vertical" data-spy = "affix" role = "navigation" aria-label = "main navigation" >
< ul >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_intro.html" > 1. Introduction< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_start.html" > 2. Getting Started< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_commands.html" > 3. Commands< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_packages.html" > 4. Packages< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_accelerate.html" > 5. Accelerating LAMMPS performance< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_howto.html" > 6. How-to discussions< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_example.html" > 7. Example problems< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_perf.html" > 8. Performance & scalability< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_tools.html" > 9. Additional tools< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_modify.html" > 10. Modifying & extending LAMMPS< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_python.html" > 11. Python interface to LAMMPS< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_errors.html" > 12. Errors< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "Section_history.html" > 13. Future and history< / a > < / li >
< / ul >
< / div >
< / nav >
< section data-toggle = "wy-nav-shift" class = "wy-nav-content-wrap" >
< nav class = "wy-nav-top" role = "navigation" aria-label = "top navigation" >
< i data-toggle = "wy-nav-top" class = "fa fa-bars" > < / i >
< a href = "Manual.html" > LAMMPS< / a >
< / nav >
< div class = "wy-nav-content" >
< div class = "rst-content" >
< div role = "navigation" aria-label = "breadcrumbs navigation" >
< ul class = "wy-breadcrumbs" >
< li > < a href = "Manual.html" > Docs< / a > » < / li >
< li > fix ti/rs command< / li >
< li class = "wy-breadcrumbs-aside" >
< a href = "http://lammps.sandia.gov" > Website< / a >
< a href = "Section_commands.html#comm" > Commands< / a >
< / li >
< / ul >
< hr / >
< / div >
< div role = "main" class = "document" itemscope = "itemscope" itemtype = "http://schema.org/Article" >
< div itemprop = "articleBody" >
< div class = "section" id = "fix-ti-rs-command" >
< span id = "index-0" > < / span > < h1 > fix ti/rs command< / h1 >
< div class = "section" id = "syntax" >
< h2 > Syntax< / h2 >
2016-09-14 23:35:03 +08:00
< pre class = "literal-block" >
fix ID group-ID ti/rs lambda_initial lambda_final t_switch t_equil keyword value ...
< / pre >
2016-05-10 01:33:12 +08:00
< ul class = "simple" >
< li > ID, group-ID are documented in < a class = "reference internal" href = "fix.html" > < span class = "doc" > fix< / span > < / a > command< / li >
< li > ti/rs = style name of this fix command< / li >
< li > lambda_initial/lambda_final = initial/final values of the coupling parameter< / li >
< li > t_switch/t_equil = number of steps of the switching/equilibration procedure< / li >
< li > keyword = < em > function< / em > < / li >
< / ul >
< pre class = "literal-block" >
< em > function< / em > value = function-ID
function-ID = ID of the switching function (1, 2 or 3)
< / pre >
< p > < strong > Example:< / strong > < / p >
< div class = "highlight-default" > < div class = "highlight" > < pre > < span > < / span > < span class = "n" > fix< / span > < span class = "n" > ref< / span > < span class = "nb" > all< / span > < span class = "n" > ti< / span > < span class = "o" > /< / span > < span class = "n" > rs< / span > < span class = "mf" > 50.0< / span > < span class = "mi" > 2000< / span > < span class = "mi" > 1000< / span >
< span class = "n" > fix< / span > < span class = "n" > vf< / span > < span class = "n" > vacancy< / span > < span class = "n" > ti< / span > < span class = "o" > /< / span > < span class = "n" > rs< / span > < span class = "mf" > 10.0< / span > < span class = "mi" > 70000< / span > < span class = "mi" > 50000< / span > < span class = "n" > function< / span > < span class = "mi" > 2< / span >
< / pre > < / div >
< / div >
< / div >
< div class = "section" id = "description" >
< h2 > Description< / h2 >
< p > This fix allows you to compute the free energy temperature dependence
by performing a thermodynamic integration procedure known as
2016-05-12 22:02:27 +08:00
Reversible Scaling < a class = "reference internal" href = "#dekoning99" > < span class = "std std-ref" > (de Koning99,< / span > < / a > < a class = "reference internal" href = "#dekoning00a" > < span class = "std std-ref" > de Koning00a)< / span > < / a > . The thermodynamic integration is performed
2016-05-10 01:33:12 +08:00
using the nonequilibrium method of Adiabatic Switching
< a class = "reference internal" href = "fix_ti_spring.html#watanabe" > < span class = "std std-ref" > (Watanabe,< / span > < / a > < a class = "reference internal" href = "fix_ti_spring.html#dekoning96" > < span class = "std std-ref" > de Koning96)< / span > < / a > .< / p >
< p > The forces on the atoms are dynamically scaled during the simulation,
the rescaling is done in the following manner:< / p >
< img alt = "_images/fix_ti_rs_force.jpg" class = "align-center" src = "_images/fix_ti_rs_force.jpg" / >
< p > where F_int is the total force on the atoms due to the interatomic
potential and lambda is the coupling parameter of the thermodynamic
integration.< / p >
< p > The fix acts as follows: during the first < em > t_equil< / em > steps after the
fix is defined the value of lambda is < em > lambda_initial< / em > , this is the
period to equilibrate the system in the lambda = < em > lambda_initial< / em >
state. After this the value of lambda changes continuously from
< em > lambda_initial< / em > to < em > lambda_final< / em > according to the function defined
using the keyword < em > function< / em > (described below), this is done in
< em > t_switch< / em > steps. Then comes the second equilibration period of
< em > t_equil< / em > to equilibrate the system in the lambda = < em > lambda_final< / em >
state. After that the switching back to the lambda = < em > lambda_initial< / em >
state is done using < em > t_switch< / em > timesteps and following the same
switching function. After this period the value of lambda is kept
equal to < em > lambda_initial< / em > indefinitely or until a < a class = "reference internal" href = "unfix.html" > < span class = "doc" > unfix< / span > < / a >
erase the fix.< / p >
< p > The description of thermodynamic integration in both directions is
done in < a class = "reference internal" href = "#dekoning00b" > < span class = "std std-ref" > de Koning00b< / span > < / a > , the main reason is to try to
eliminate the dissipated heat due to the nonequilibrium process.< / p >
< p > The < em > function< / em > keyword allows the use of three different switching
rates. The option < em > 1< / em > results in a constant rescaling where the lambda
parameter changes at a constant rate during the switching time
according to the switching function< / p >
< img alt = "_images/fix_ti_rs_function_1.jpg" class = "align-center" src = "_images/fix_ti_rs_function_1.jpg" / >
< p > where tau is the scaled time variable t/t_switch. This switching
function has the characteristic that the temperature scaling is faster
at temperatures closer to the final temperature of the procedure. The
option number < em > 2< / em > performs the switching at a rate defined by the
following switching function< / p >
< img alt = "_images/fix_ti_rs_function_2.jpg" class = "align-center" src = "_images/fix_ti_rs_function_2.jpg" / >
< p > This switching function has the characteristic that the temperature
scaling occurs at a constant rate during all the procedure. The option
number < em > 3< / em > performs the switching at a rate defined by the following
switching function< / p >
< img alt = "_images/fix_ti_rs_function_3.jpg" class = "align-center" src = "_images/fix_ti_rs_function_3.jpg" / >
< p > This switching function has the characteristic that the temperature
scaling is faster at temperatures closer to the initial temperature of
the procedure.< / p >
< p > An example script using this command is provided in the
examples/USER/misc/ti directory.< / p >
2016-09-14 23:35:03 +08:00
< p > < strong > Restart, fix_modify, output, run start/stop, minimize info:< / strong > < / p >
2016-05-10 01:33:12 +08:00
< p > No information about this fix is written to < a class = "reference internal" href = "restart.html" > < span class = "doc" > binary restart files< / span > < / a > .< / p >
< p > This fix computes a global vector quantitie which can be accessed by
various < a class = "reference internal" href = "Section_howto.html#howto-15" > < span class = "std std-ref" > output commands< / span > < / a > . The vector has
2 positions, the first one is the coupling parameter lambda and the
second one is the time derivative of lambda. The scalar and vector
values calculated by this fix are “ extensive” .< / p >
< p > No parameter of this fix can be used with the < em > start/stop< / em > keywords of
the < a class = "reference internal" href = "run.html" > < span class = "doc" > run< / span > < / a > command.< / p >
< p > The forces due to this fix are imposed during an energy minimization,
invoked by the < a class = "reference internal" href = "minimize.html" > < span class = "doc" > minimize< / span > < / a > command.< / p >
< / div >
< div class = "section" id = "related-commands" >
< h2 > Related commands< / h2 >
< p > < a class = "reference internal" href = "fix_ti_spring.html" > < span class = "doc" > fix ti/spring< / span > < / a > < / p >
< / div >
< div class = "section" id = "restrictions" >
< h2 > Restrictions< / h2 >
< p > This command is part of the USER-MISC package. It is only enabled if
LAMMPS was built with those packages. See the < a class = "reference internal" href = "Section_start.html#start-3" > < span class = "std std-ref" > Making LAMMPS< / span > < / a > section for more info.< / p >
< / div >
< div class = "section" id = "default" >
< h2 > Default< / h2 >
< p > The keyword default is function = 1.< / p >
< hr class = "docutils" / >
< p id = "dekoning99" > < strong > (de Koning 99)< / strong > M. de Koning, A. Antonelli and S. Yip, Phys Rev Lett, 83, 3973 (1999).< / p >
< p id = "watanabe" > < strong > (Watanabe)< / strong > M. Watanabe and W. P. Reinhardt, Phys Rev Lett, 65, 3301 (1990).< / p >
< p id = "dekoning96" > < strong > (de Koning 96)< / strong > M. de Koning and A. Antonelli, Phys Rev E, 53, 465 (1996).< / p >
< p id = "dekoning00a" > < strong > (de Koning 00a)< / strong > M. de Koning, A. Antonelli and S. Yip, J Chem Phys, 115, 11025 (2000).< / p >
< p id = "dekoning00b" > < strong > (de Koning 00b)< / strong > M. de Koning et al., Computing in Science & Engineering, 2, 88 (2000).< / p >
< / div >
< / div >
< / div >
< / div >
< footer >
< hr / >
< div role = "contentinfo" >
< p >
© Copyright 2013 Sandia Corporation.
< / p >
< / div >
Built with < a href = "http://sphinx-doc.org/" > Sphinx< / a > using a < a href = "https://github.com/snide/sphinx_rtd_theme" > theme< / a > provided by < a href = "https://readthedocs.org" > Read the Docs< / a > .
< / footer >
< / div >
< / div >
< / section >
< / div >
< script type = "text/javascript" >
var DOCUMENTATION_OPTIONS = {
URL_ROOT:'./',
VERSION:'',
COLLAPSE_INDEX:false,
FILE_SUFFIX:'.html',
HAS_SOURCE: true
};
< / script >
< script type = "text/javascript" src = "_static/jquery.js" > < / script >
< script type = "text/javascript" src = "_static/underscore.js" > < / script >
< script type = "text/javascript" src = "_static/doctools.js" > < / script >
< script type = "text/javascript" src = "https://cdn.mathjax.org/mathjax/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML" > < / script >
< script type = "text/javascript" src = "_static/sphinxcontrib-images/LightBox2/lightbox2/js/jquery-1.11.0.min.js" > < / script >
< script type = "text/javascript" src = "_static/sphinxcontrib-images/LightBox2/lightbox2/js/lightbox.min.js" > < / script >
< script type = "text/javascript" src = "_static/sphinxcontrib-images/LightBox2/lightbox2-customize/jquery-noconflict.js" > < / script >
< script type = "text/javascript" src = "_static/js/theme.js" > < / script >
< script type = "text/javascript" >
jQuery(function () {
SphinxRtdTheme.StickyNav.enable();
});
< / script >
< / body >
< / html >