lammps/doc/fix_temp_berendsen.html

345 lines
17 KiB
HTML

<!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 temp/berendsen command &mdash; 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 &amp; 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 &amp; 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>
&nbsp;
</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> &raquo;</li>
<li>fix temp/berendsen 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-temp-berendsen-command">
<span id="index-0"></span><h1>fix temp/berendsen command<a class="headerlink" href="#fix-temp-berendsen-command" title="Permalink to this headline"></a></h1>
</div>
<div class="section" id="fix-temp-berendsen-cuda-command">
<h1>fix temp/berendsen/cuda command<a class="headerlink" href="#fix-temp-berendsen-cuda-command" title="Permalink to this headline"></a></h1>
<div class="section" id="syntax">
<h2>Syntax<a class="headerlink" href="#syntax" title="Permalink to this headline"></a></h2>
<div class="highlight-python"><div class="highlight"><pre>fix ID group-ID temp/berendsen Tstart Tstop Tdamp
</pre></div>
</div>
<ul class="simple">
<li>ID, group-ID are documented in <a class="reference internal" href="fix.html"><em>fix</em></a> command</li>
<li>temp/berendsen = style name of this fix command</li>
<li>Tstart,Tstop = desired temperature at start/end of run</li>
</ul>
<div class="highlight-python"><div class="highlight"><pre>Tstart can be a variable (see below)
</pre></div>
</div>
<ul class="simple">
<li>Tdamp = temperature damping parameter (time units)</li>
</ul>
</div>
<div class="section" id="examples">
<h2>Examples<a class="headerlink" href="#examples" title="Permalink to this headline"></a></h2>
<div class="highlight-python"><div class="highlight"><pre>fix 1 all temp/berendsen 300.0 300.0 100.0
</pre></div>
</div>
</div>
<div class="section" id="description">
<h2>Description<a class="headerlink" href="#description" title="Permalink to this headline"></a></h2>
<p>Reset the temperature of a group of atoms by using a Berendsen
thermostat <a class="reference internal" href="#berendsen"><span>(Berendsen)</span></a>, which rescales their velocities
every timestep.</p>
<p>The thermostat is applied to only the translational degrees of freedom
for the particles, which is an important consideration for finite-size
particles which have rotational degrees of freedom are being
thermostatted with this fix. The translational degrees of freedom can
also have a bias velocity removed from them before thermostatting
takes place; see the description below.</p>
<p>The desired temperature at each timestep is a ramped value during the
run from <em>Tstart</em> to <em>Tstop</em>. The <em>Tdamp</em> parameter is specified in
time units and determines how rapidly the temperature is relaxed. For
example, a value of 100.0 means to relax the temperature in a timespan
of (roughly) 100 time units (tau or fmsec or psec - see the
<a class="reference internal" href="units.html"><em>units</em></a> command).</p>
<p><em>Tstart</em> can be specified as an equal-style <a class="reference internal" href="variable.html"><em>variable</em></a>.
In this case, the <em>Tstop</em> setting is ignored. If the value is a
variable, it should be specified as v_name, where name is the variable
name. In this case, the variable will be evaluated each timestep, and
its value used to determine the target temperature.</p>
<p>Equal-style variables can specify formulas with various mathematical
functions, and include <a class="reference internal" href="thermo_style.html"><em>thermo_style</em></a> command
keywords for the simulation box parameters and timestep and elapsed
time. Thus it is easy to specify a time-dependent temperature.</p>
<div class="admonition note">
<p class="first admonition-title">Note</p>
<p class="last">Unlike the <a class="reference internal" href="fix_nh.html"><em>fix nvt</em></a> command which performs
Nose/Hoover thermostatting AND time integration, this fix does NOT
perform time integration. It only modifies velocities to effect
thermostatting. Thus you must use a separate time integration fix,
like <a class="reference internal" href="fix_nve.html"><em>fix nve</em></a> to actually update the positions of atoms
using the modified velocities. Likewise, this fix should not normally
be used on atoms that also have their temperature controlled by
another fix - e.g. by <a class="reference internal" href="fix_nh.html"><em>fix nvt</em></a> or <a class="reference internal" href="fix_langevin.html"><em>fix langevin</em></a> commands.</p>
</div>
<p>See <a class="reference internal" href="Section_howto.html#howto-16"><span>this howto section</span></a> of the manual for
a discussion of different ways to compute temperature and perform
thermostatting.</p>
<p>This fix computes a temperature each timestep. To do this, the fix
creates its own compute of style &#8220;temp&#8221;, as if this command had been
issued:</p>
<div class="highlight-python"><div class="highlight"><pre>compute fix-ID_temp group-ID temp
</pre></div>
</div>
<p>See the <a class="reference internal" href="compute_temp.html"><em>compute temp</em></a> command for details. Note
that the ID of the new compute is the fix-ID + underscore + &#8220;temp&#8221;,
and the group for the new compute is the same as the fix group.</p>
<p>Note that this is NOT the compute used by thermodynamic output (see
the <a class="reference internal" href="thermo_style.html"><em>thermo_style</em></a> command) with ID = <em>thermo_temp</em>.
This means you can change the attributes of this fix&#8217;s temperature
(e.g. its degrees-of-freedom) via the
<a class="reference internal" href="compute_modify.html"><em>compute_modify</em></a> command or print this temperature
during thermodynamic output via the <a class="reference internal" href="thermo_style.html"><em>thermo_style custom</em></a> command using the appropriate compute-ID.
It also means that changing attributes of <em>thermo_temp</em> will have no
effect on this fix.</p>
<p>Like other fixes that perform thermostatting, this fix can be used
with <a class="reference internal" href="compute.html"><em>compute commands</em></a> that calculate a temperature
after removing a &#8220;bias&#8221; from the atom velocities. E.g. removing the
center-of-mass velocity from a group of atoms or only calculating
temperature on the x-component of velocity or only calculating
temperature for atoms in a geometric region. This is not done by
default, but only if the <a class="reference internal" href="fix_modify.html"><em>fix_modify</em></a> command is used
to assign a temperature compute to this fix that includes such a bias
term. See the doc pages for individual <a class="reference internal" href="compute.html"><em>compute commands</em></a> to determine which ones include a bias. In
this case, the thermostat works in the following manner: the current
temperature is calculated taking the bias into account, bias is
removed from each atom, thermostatting is performed on the remaining
thermal degrees of freedom, and the bias is added back in.</p>
<hr class="docutils" />
<p>Styles with a <em>cuda</em> suffix are functionally the same as the
corresponding style without the suffix. They have been optimized to
run faster, depending on your available hardware, as discussed in
<a class="reference internal" href="Section_accelerate.html"><em>Section_accelerate</em></a> of the manual. The
accelerated styles take the same arguments and should produce the same
results, except for round-off and precision issues.</p>
<p>These accelerated styles are part of the USER-CUDA package. They are
only enabled if LAMMPS was built with that package. See the <a class="reference internal" href="Section_start.html#start-3"><span>Making LAMMPS</span></a> section for more info.</p>
<p>You can specify the accelerated styles explicitly in your input script
by including their suffix, or you can use the <a class="reference internal" href="Section_start.html#start-7"><span>-suffix command-line switch</span></a> when you invoke LAMMPS, or you can
use the <a class="reference internal" href="suffix.html"><em>suffix</em></a> command in your input script.</p>
<p>See <a class="reference internal" href="Section_accelerate.html"><em>Section_accelerate</em></a> of the manual for
more instructions on how to use the accelerated styles effectively.</p>
</div>
<hr class="docutils" />
<div class="section" id="restart-fix-modify-output-run-start-stop-minimize-info">
<h2>Restart, fix_modify, output, run start/stop, minimize info<a class="headerlink" href="#restart-fix-modify-output-run-start-stop-minimize-info" title="Permalink to this headline"></a></h2>
<p>No information about this fix is written to <a class="reference internal" href="restart.html"><em>binary restart files</em></a>.</p>
<p>The <a class="reference internal" href="fix_modify.html"><em>fix_modify</em></a> <em>temp</em> option is supported by this
fix. You can use it to assign a temperature <a class="reference internal" href="compute.html"><em>compute</em></a>
you have defined to this fix which will be used in its thermostatting
procedure, as described above. For consistency, the group used by
this fix and by the compute should be the same.</p>
<p>The <a class="reference internal" href="fix_modify.html"><em>fix_modify</em></a> <em>energy</em> option is supported by this
fix to add the energy change implied by a velocity rescaling to the
system&#8217;s potential energy as part of <a class="reference internal" href="thermo_style.html"><em>thermodynamic output</em></a>.</p>
<p>This fix computes a global scalar which can be accessed by various
<a class="reference internal" href="Section_howto.html#howto-15"><span>output commands</span></a>. The scalar is the
cummulative energy change due to this fix. The scalar value
calculated by this fix is &#8220;extensive&#8221;.</p>
<p>This fix can ramp its target temperature over multiple runs, using the
<em>start</em> and <em>stop</em> keywords of the <a class="reference internal" href="run.html"><em>run</em></a> command. See the
<a class="reference internal" href="run.html"><em>run</em></a> command for details of how to do this.</p>
<p>This fix is not invoked during <a class="reference internal" href="minimize.html"><em>energy minimization</em></a>.</p>
</div>
<div class="section" id="restrictions">
<h2>Restrictions<a class="headerlink" href="#restrictions" title="Permalink to this headline"></a></h2>
<p>This fix can be used with dynamic groups as defined by the
<a class="reference internal" href="group.html"><em>group</em></a> command. Likewise it can be used with groups to
which atoms are added or deleted over time, e.g. a deposition
simulation. However, the conservation properties of the thermostat
and barostat are defined for systems with a static set of atoms. You
may observe odd behavior if the atoms in a group vary dramatically
over time or the atom count becomes very small.</p>
</div>
<div class="section" id="related-commands">
<h2>Related commands<a class="headerlink" href="#related-commands" title="Permalink to this headline"></a></h2>
<p><a class="reference internal" href="fix_nve.html"><em>fix nve</em></a>, <a class="reference internal" href="fix_nh.html"><em>fix nvt</em></a>, <a class="reference internal" href="fix_temp_rescale.html"><em>fix temp/rescale</em></a>, <a class="reference internal" href="fix_langevin.html"><em>fix langevin</em></a>,
<a class="reference internal" href="fix_modify.html"><em>fix_modify</em></a>, <a class="reference internal" href="compute_temp.html"><em>compute temp</em></a>,
<a class="reference internal" href="fix_press_berendsen.html"><em>fix press/berendsen</em></a></p>
<p><strong>Default:</strong> none</p>
<hr class="docutils" />
<p id="berendsen"><strong>(Berendsen)</strong> Berendsen, Postma, van Gunsteren, DiNola, Haak, J Chem
Phys, 81, 3684 (1984).</p>
</div>
</div>
</div>
</div>
<footer>
<hr/>
<div role="contentinfo">
<p>
&copy; 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>