Page MenuHomec4science

fix_rx.html
No OneTemporary

File Metadata

Created
Sat, Aug 10, 18:19

fix_rx.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 rx 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 rx 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-rx-command">
<span id="index-0"></span><h1>fix rx command</h1>
<div class="section" id="syntax">
<h2>Syntax</h2>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">fix</span> <span class="n">ID</span> <span class="n">group</span><span class="o">-</span><span class="n">ID</span> <span class="n">rx</span> <span class="n">file</span> <span class="n">localTemp</span> <span class="n">matrix</span> <span class="n">solver</span> <span class="n">minSteps</span> <span class="o">...</span>
</pre></div>
</div>
<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>rx = style name of this fix command</li>
<li>file = filename containing the reaction kinetic equations and Arrhenius parameters</li>
<li>localTemp = <em>none,lucy</em> = no local temperature averaging or local temperature defined through Lucy weighting function</li>
<li>matrix = <em>sparse, dense</em> format for the stoichiometric matrix</li>
<li>solver = <em>lammps_rk4,rkf45</em> = rk4 is an explicit 4th order Runge-Kutta method; rkf45 is an adaptive 4th-order Runge-Kutta-Fehlberg method</li>
<li>minSteps = # of steps for rk4 solver or minimum # of steps for rkf45 (rk4 or rkf45)</li>
<li>maxSteps = maximum number of steps for the rkf45 solver (rkf45 only)</li>
<li>relTol = relative tolerance for the rkf45 solver (rkf45 only)</li>
<li>absTol = absolute tolernace for the rkf45 solver (rkf45 only)</li>
<li>diag = Diagnostics frequency for the rkf45 solver (optional, rkf45 only)</li>
</ul>
</div>
<div class="section" id="examples">
<h2>Examples</h2>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">fix</span> <span class="mi">1</span> <span class="nb">all</span> <span class="n">rx</span> <span class="n">kinetics</span><span class="o">.</span><span class="n">rx</span> <span class="n">none</span> <span class="n">dense</span> <span class="n">lammps_rk4</span>
<span class="n">fix</span> <span class="mi">1</span> <span class="nb">all</span> <span class="n">rx</span> <span class="n">kinetics</span><span class="o">.</span><span class="n">rx</span> <span class="n">none</span> <span class="n">sparse</span> <span class="n">lammps_rk4</span> <span class="mi">1</span>
<span class="n">fix</span> <span class="mi">1</span> <span class="nb">all</span> <span class="n">rx</span> <span class="n">kinetics</span><span class="o">.</span><span class="n">rx</span> <span class="n">lucy</span> <span class="n">sparse</span> <span class="n">lammps_rk4</span> <span class="mi">10</span>
<span class="n">fix</span> <span class="mi">1</span> <span class="nb">all</span> <span class="n">rx</span> <span class="n">kinetics</span><span class="o">.</span><span class="n">rx</span> <span class="n">none</span> <span class="n">dense</span> <span class="n">rkf45</span> <span class="mi">1</span> <span class="mi">100</span> <span class="mi">1</span><span class="n">e</span><span class="o">-</span><span class="mi">6</span> <span class="mi">1</span><span class="n">e</span><span class="o">-</span><span class="mi">8</span>
<span class="n">fix</span> <span class="mi">1</span> <span class="nb">all</span> <span class="n">rx</span> <span class="n">kinetics</span><span class="o">.</span><span class="n">rx</span> <span class="n">none</span> <span class="n">dense</span> <span class="n">rkf45</span> <span class="mi">1</span> <span class="mi">100</span> <span class="mi">1</span><span class="n">e</span><span class="o">-</span><span class="mi">6</span> <span class="mi">1</span><span class="n">e</span><span class="o">-</span><span class="mi">8</span> <span class="o">-</span><span class="mi">1</span>
</pre></div>
</div>
</div>
<div class="section" id="description">
<h2>Description</h2>
<p>Fix <em>rx</em> solves the reaction kinetic ODEs for a given reaction set that is
defined within the file associated with this command.</p>
<p>For a general reaction such that</p>
<img alt="_images/fix_rx_reaction.jpg" class="align-center" src="_images/fix_rx_reaction.jpg" />
<p>the reaction rate equation is defined to be of the form</p>
<img alt="_images/fix_rx_reactionRate.jpg" class="align-center" src="_images/fix_rx_reactionRate.jpg" />
<p>In the current implementation, the exponents are defined to be equal
to the stoichiometric coefficients. A given reaction set consisting
of <em>n</em> reaction equations will contain a total of <em>m</em> species. A set
of <em>m</em> ordinary differential equations (ODEs) that describe the change
in concentration of a given species as a function of time are then
constructed based on the <em>n</em> reaction rate equations.</p>
<p>The ODE systems are solved over the full DPD timestep <em>dt</em> using either a 4th
order Runge-Kutta <em>rk4</em> method with a fixed step-size <em>h</em>, specified
by the <em>lammps_rk4</em> keyword, or a 4th order Runge-Kutta-Fehlberg (rkf45) method
with an adaptive step-size for <em>h</em>. The number of ODE steps per DPD timestep
for the rk4 method is optionally specified immediately after the rk4
keyword. The ODE step-size is set as <em>dt/num_steps</em>. Smaller
step-sizes tend to yield more accurate results but there is not
control on the error. For error control, use the rkf45 ODE solver.</p>
<p>The rkf45 method adjusts the step-size so that the local truncation error is held
within the specified absolute and relative tolerances. The initial step-size <em>h0</em>
can be specified by the user or estimated internally. It is recommeded that the user
specify <em>h0</em> since this will generally reduced the number of ODE integration steps
required. <em>h0</em> is defined as <em>dt / min_steps</em> if min_steps &gt;= 1. If min_steps == 0,
<em>h0</em> is estimated such that an explicit Euler method would likely produce
an acceptable solution. This is generally overly conservative for the 4th-order
method and users are advised to specify <em>h0</em> as some fraction of the DPD timestep.
For small DPD timesteps, only one step may be necessary depending upon the tolerances.
Note that more than min_steps ODE steps may be taken depending upon the ODE stiffness
but no more than max_steps will be taken. If max_steps is reached, an error warning
is printed and the simulation is stopped.</p>
<p>After each ODE step, the solution error <em>e</em> is tested and weighted using the absTol
and relTol values. The error vector is weighted as <em>e</em> / (relTol * <a href="#id1"><span class="problematic" id="id2">|</span></a><em>u</em>| + absTol)
where <em>u</em> is the solution vector. If the norm of the error is &lt;= 1, the solution is
accepted, <em>h</em> is increased by a proportional amount, and the next ODE step is begun.
Otherwise, <em>h</em> is shrunk and the ODE step is repeated.</p>
<p>Run-time diagnostics are available for the rkf45 ODE solver. The frequency
(in time-steps) that diagnostics are reported is controlled by the last (optional)
12th argument. A negative frequency means that diagnostics are reported once at the
end of each run. A positive value N means that the diagnostics are reported once
per N time-steps.</p>
<p>The diagnostics report the average # of integrator steps and RHS function evaluations
and run-time per ODE as well as the the average/RMS/min/max per process. If the
reporting frequency is 1, the RMS/min/max per ODE are also reported. The per ODE
statistics can be used to adjust the tolerance and min/max step parameters. The
statistics per MPI process can be useful to examine any load imbalance caused by the
adaptive ODE solver. (Some DPD particles can take longer to solve than others. This
can lead to an imbalance across the MPI processes.)</p>
<hr class="docutils" />
<p>The filename specifies a file that contains the entire set of reaction
kinetic equations and corresponding Arrhenius parameters. The format of
this file is described below.</p>
<p>There is no restriction on the total number or reaction equations that
are specified. The species names are arbitrary string names that are
associated with the species concentrations. Each species in a given
reaction must be preceded by it&#8217;s stoichiometric coefficient. The
only delimiters that are recognized between the species are either a
<em>+</em> or <em>=</em> character. The <em>=</em> character corresponds to an
irreversible reaction. After specifying the reaction, the reaction
rate constant is determined through the temperature dependent
Arrhenius equation:</p>
<img alt="_images/fix_rx.jpg" class="align-center" src="_images/fix_rx.jpg" />
<p>where <em>A</em> is the Arrhenius factor in time units or concentration/time
units, <em>n</em> is the unitless exponent of the temperature dependence, and
<em>E_a</em> is the activation energy in energy units. The temperature
dependence can be removed by specifying the exponent as zero.</p>
<p>The internal temperature of the coarse-grained particles can be used
in constructing the reaction rate constants at every DPD timestep by
specifying the keyword <em>none</em>. Alternatively, the keyword <em>lucy</em> can
be specified to compute a local-average particle internal temperature
for use in the reaction rate constant expressions. The local-average
particle internal temperature is defined as:</p>
<img alt="_images/fix_rx_localTemp.jpg" class="align-center" src="_images/fix_rx_localTemp.jpg" />
<p>where the Lucy function is expressed as:</p>
<img alt="_images/fix_rx_localTemp2.jpg" class="align-center" src="_images/fix_rx_localTemp2.jpg" />
<p>The self-particle interaction is included in the above equation.</p>
<p>The stoichiometric coefficients for the reaction mechanism are stored
in either a sparse or dense matrix format. The dense matrix should only be
used for small reaction mechanisms. The sparse matrix should be used when there
are many reactions (e.g., more than 5). This allows the number of reactions and
species to grow while keeping the computational cost tractable. The matrix
format can be specified as using either the <em>sparse</em> or <em>dense</em> keywords.
If all stoichiometric coefficients for a reaction are small integers (whole
numbers &lt;= 3), a fast exponential function is used. This can save significant
computational time so users are encouraged to use integer coefficients
where possible.</p>
<hr class="docutils" />
<p>The format of a tabulated file is as follows (without the
parenthesized comments):</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="c1"># Rxn equations and parameters (one or more comment or blank lines)</span>
</pre></div>
</div>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="mf">1.0</span> <span class="n">hcn</span> <span class="o">+</span> <span class="mf">1.0</span> <span class="n">no2</span> <span class="o">=</span> <span class="mf">1.0</span> <span class="n">no</span> <span class="o">+</span> <span class="mf">0.5</span> <span class="n">n2</span> <span class="o">+</span> <span class="mf">0.5</span> <span class="n">h2</span> <span class="o">+</span> <span class="mf">1.0</span> <span class="n">co</span> <span class="mf">2.49E+01</span> <span class="mf">0.0</span> <span class="mf">1.34</span> <span class="p">(</span><span class="n">rxn</span> <span class="n">equation</span><span class="p">,</span> <span class="n">A</span><span class="p">,</span> <span class="n">n</span><span class="p">,</span> <span class="n">Ea</span><span class="p">)</span>
<span class="mf">1.0</span> <span class="n">hcn</span> <span class="o">+</span> <span class="mf">1.0</span> <span class="n">no</span> <span class="o">=</span> <span class="mf">1.0</span> <span class="n">co</span> <span class="o">+</span> <span class="mf">1.0</span> <span class="n">n2</span> <span class="o">+</span> <span class="mf">0.5</span> <span class="n">h2</span> <span class="mf">2.16E+00</span> <span class="mf">0.0</span> <span class="mf">1.52</span>
<span class="o">...</span>
<span class="mf">1.0</span> <span class="n">no</span> <span class="o">+</span> <span class="mf">1.0</span> <span class="n">co</span> <span class="o">=</span> <span class="mf">0.5</span> <span class="n">n2</span> <span class="o">+</span> <span class="mf">1.0</span> <span class="n">co2</span> <span class="mf">1.66E+06</span> <span class="mf">0.0</span> <span class="mf">0.69</span>
</pre></div>
</div>
<p>A section begins with a non-blank line whose 1st character is not a
&#8220;#&#8221;; blank lines or lines starting with &#8220;#&#8221; can be used as comments
between sections.</p>
<p>Following a blank line, the next N lines list the N reaction
equations. Each species within the reaction equation is specified
through its stoichiometric coefficient and a species tag. Reactant
species are specified on the left-hand side of the equation and
product species are specified on the right-hand side of the equation.
After specifying the reactant and product species, the final three
arguments of each line represent the Arrhenius parameter <em>A</em>, the
temperature exponent <em>n</em>, and the activation energy <em>Ea</em>.</p>
<p>Note that the species tags that are defined in the reaction equations
are used by the <a class="reference internal" href="fix_eos_table_rx.html"><span class="doc">fix eos/table/rx</span></a> command to
define the thermodynamic properties of each species. Furthermore, the
number of species molecules (i.e., concentration) can be specified
either with the <a class="reference internal" href="set.html"><span class="doc">set</span></a> command using the &#8220;<a href="#id3"><span class="problematic" id="id4">d_</span></a>&#8221; prefix or by
reading directly the concentrations from a data file. For the latter
case, the <a class="reference internal" href="read_data.html"><span class="doc">read_data</span></a> command with the fix keyword
should be specified, where the fix-ID will be the &#8220;fix rx`ID with a &lt;SPECIES&#8221;&gt;`_ suffix, e.g.</p>
<p>fix foo all rx reaction.file ...
read_data data.dpd fix foo_SPECIES NULL Species</p>
</div>
<hr class="docutils" />
<div class="section" id="restrictions">
<h2>Restrictions</h2>
<p>This command is part of the USER-DPD package. It is only enabled if
LAMMPS was built with that package. 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>
<p>This command also requires use of the <a class="reference internal" href="atom_style.html"><span class="doc">atom_style dpd</span></a>
command.</p>
<p>This command can only be used with a constant energy or constant
enthalpy DPD simulation.</p>
</div>
<div class="section" id="related-commands">
<h2>Related commands</h2>
<p><a class="reference internal" href="fix_eos_table_rx.html"><span class="doc">fix eos/table/rx</span></a>,
<a class="reference internal" href="fix_shardlow.html"><span class="doc">fix shardlow</span></a>,
<span class="xref doc">pair dpd/fdt/energy</span></p>
<p><strong>Default:</strong> none</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>

Event Timeline