lammps/doc/html/body.html

452 lines
24 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>Body particles &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>Body particles</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="body-particles">
<h1>Body particles</h1>
<p><strong>Overview:</strong></p>
<p>This doc page is not about a LAMMPS input script command, but about
body particles, which are generalized finite-size particles.
Individual body particles can represent complex entities, such as
surface meshes of discrete points, collections of sub-particles,
deformable objects, etc. Note that other kinds of finite-size
spherical and aspherical particles are also supported by LAMMPS, such
as spheres, ellipsoids, line segments, and triangles, but they are
simpler entities that body particles. See <a class="reference internal" href="Section_howto.html#howto-14"><span class="std std-ref">Section_howto 14</span></a> for a general overview of all these
particle types.</p>
<p>Body particles are used via the <a class="reference internal" href="atom_style.html"><span class="doc">atom_style body</span></a>
command. It takes a body style as an argument. The current body
styles supported by LAMMPS are as follows. The name in the first
column is used as the <em>bstyle</em> argument for the <a class="reference internal" href="atom_style.html"><span class="doc">atom_style body</span></a> command.</p>
<table border="1" class="docutils">
<colgroup>
<col width="35%" />
<col width="65%" />
</colgroup>
<tbody valign="top">
<tr class="row-odd"><td><em>nparticle</em></td>
<td>rigid body with N sub-particles</td>
</tr>
<tr class="row-even"><td><em>rounded/polygon</em></td>
<td>2d convex polygon with N vertices</td>
</tr>
</tbody>
</table>
<p>The body style determines what attributes are stored for each body and
thus how they can be used to compute pairwise body/body or
bond/non-body (point particle) interactions. More details of each
style are described below.</p>
<div class="admonition note">
<p class="first admonition-title">Note</p>
<p class="last">The rounded/polygon style listed in the table above and
described below has not yet been relesed in LAMMPS. It will be soon.</p>
</div>
<p>We hope to add more styles in the future. See <a class="reference internal" href="Section_modify.html#mod-12"><span class="std std-ref">Section_modify 12</span></a> for details on how to add a new body
style to the code.</p>
<hr class="docutils" />
<p><strong>When to use body particles:</strong></p>
<p>You should not use body particles to model a rigid body made of
simpler particles (e.g. point, sphere, ellipsoid, line segment,
triangular particles), if the interaction between pairs of rigid
bodies is just the summation of pairwise interactions between the
simpler particles. LAMMPS already supports this kind of model via the
<a class="reference internal" href="fix_rigid.html"><span class="doc">fix rigid</span></a> command. Any of the numerous pair styles
that compute interactions between simpler particles can be used. The
<a class="reference internal" href="fix_rigid.html"><span class="doc">fix rigid</span></a> command time integrates the motion of the
rigid bodies. All of the standard LAMMPS commands for thermostatting,
adding constraints, performing output, etc will operate as expected on
the simple particles.</p>
<p>By contrast, when body particles are used, LAMMPS treats an entire
body as a single particle for purposes of computing pairwise
interactions, building neighbor lists, migrating particles between
processors, outputting particles to a dump file, etc. This means that
interactions between pairs of bodies or between a body and non-body
(point) particle need to be encoded in an appropriate pair style. If
such a pair style were to mimic the <a class="reference internal" href="fix_rigid.html"><span class="doc">fix rigid</span></a> model,
it would need to loop over the entire collection of interactions
between pairs of simple particles within the two bodies, each time a
single body/body interaction was computed.</p>
<p>Thus it only makes sense to use body particles and develop such a pair
style, when particle/particle interactions are more complex than what
the <a class="reference internal" href="fix_rigid.html"><span class="doc">fix rigid</span></a> command can already calculate. For
example, if particles have one or more of the following attributes:</p>
<ul class="simple">
<li>represented by a surface mesh</li>
<li>represented by a collection of geometric entities (e.g. planes + spheres)</li>
<li>deformable</li>
<li>internal stress that induces fragmentation</li>
</ul>
<p>then the interaction between pairs of particles is likely to be more
complex than the summation of simple sub-particle interactions. An
example is contact or frictional forces between particles with planar
sufaces that inter-penetrate.</p>
<p>These are additional LAMMPS commands that can be used with body
particles of different styles</p>
<table border="1" class="docutils">
<colgroup>
<col width="48%" />
<col width="52%" />
</colgroup>
<tbody valign="top">
<tr class="row-odd"><td><a class="reference internal" href="fix_nve_body.html"><span class="doc">fix nve/body</span></a></td>
<td>integrate motion of a body particle in NVE ensemble</td>
</tr>
<tr class="row-even"><td><a class="reference internal" href="fix_nvt_body.html"><span class="doc">fix nvt/body</span></a></td>
<td>ditto for NVT ensemble</td>
</tr>
<tr class="row-odd"><td><a class="reference internal" href="fix_npt_body.html"><span class="doc">fix npt/body</span></a></td>
<td>ditto for NPT ensemble</td>
</tr>
<tr class="row-even"><td><a class="reference internal" href="fix_nph_body.html"><span class="doc">fix nph/body</span></a></td>
<td>ditto for NPH ensemble</td>
</tr>
<tr class="row-odd"><td><a class="reference internal" href="compute_body_local.html"><span class="doc">compute body/local</span></a></td>
<td>store sub-particle attributes of a body particle</td>
</tr>
<tr class="row-even"><td><a class="reference internal" href="compute_temp_body.html"><span class="doc">compute temp/body</span></a></td>
<td>compute temperature of body particles</td>
</tr>
<tr class="row-odd"><td><a class="reference internal" href="dump.html"><span class="doc">dump local</span></a></td>
<td>output sub-particle attributes of a body particle</td>
</tr>
<tr class="row-even"><td><a class="reference internal" href="dump_image.html"><span class="doc">dump image</span></a></td>
<td>output body particle attributes as an image</td>
</tr>
</tbody>
</table>
<p>The pair styles defined for use with specific body styles are listed
in the sections below.</p>
<hr class="docutils" />
<p><strong>Specifics of body style nparticle:</strong></p>
<p>The <em>nparticle</em> body style represents body particles as a rigid body
with a variable number N of sub-particles. It is provided as a
vanillia, prototypical example of a body particle, although as
mentioned above, the <a class="reference internal" href="fix_rigid.html"><span class="doc">fix rigid</span></a> command already
duplicates its functionality.</p>
<p>The atom_style body command for this body style takes two additional
arguments:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">atom_style</span> <span class="n">body</span> <span class="n">nparticle</span> <span class="n">Nmin</span> <span class="n">Nmax</span>
<span class="n">Nmin</span> <span class="o">=</span> <span class="n">minimum</span> <span class="c1"># of sub-particles in any body in the system</span>
<span class="n">Nmax</span> <span class="o">=</span> <span class="n">maximum</span> <span class="c1"># of sub-particles in any body in the system</span>
</pre></div>
</div>
<p>The Nmin and Nmax arguments are used to bound the size of data
structures used internally by each particle.</p>
<p>When the <a class="reference internal" href="read_data.html"><span class="doc">read_data</span></a> command reads a data file for this
body style, the following information must be provided for each entry
in the <em>Bodies</em> section of the data file:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">atom</span><span class="o">-</span><span class="n">ID</span> <span class="mi">1</span> <span class="n">M</span>
<span class="n">N</span>
<span class="n">ixx</span> <span class="n">iyy</span> <span class="n">izz</span> <span class="n">ixy</span> <span class="n">ixz</span> <span class="n">iyz</span>
<span class="n">x1</span> <span class="n">y1</span> <span class="n">z1</span>
<span class="o">...</span>
<span class="n">xN</span> <span class="n">yN</span> <span class="n">zN</span>
</pre></div>
</div>
<p>N is the number of sub-particles in the body particle. M = 6 + 3*N.
The integer line has a single value N. The floating point line(s)
list 6 moments of inertia followed by the coordinates of the N
sub-particles (x1 to zN) as 3N values. These values can be listed on
as many lines as you wish; see the <a class="reference internal" href="read_data.html"><span class="doc">read_data</span></a> command
for more details.</p>
<p>The 6 moments of inertia (ixx,iyy,izz,ixy,ixz,iyz) should be the
values consistent with the current orientation of the rigid body
around its center of mass. The values are with respect to the
simulation box XYZ axes, not with respect to the prinicpal axes of the
rigid body itself. LAMMPS performs the latter calculation internally.
The coordinates of each sub-particle are specified as its x,y,z
displacement from the center-of-mass of the body particle. The
center-of-mass position of the particle is specified by the x,y,z
values in the <em>Atoms</em> section of the data file, as is the total mass
of the body particle.</p>
<p>The <a class="reference internal" href="pair_body.html"><span class="doc">pair_style body</span></a> command can be used with this
body style to compute body/body and body/non-body interactions.</p>
<p>For output purposes via the <a class="reference internal" href="compute_body_local.html"><span class="doc">compute body/local</span></a> and <a class="reference internal" href="dump.html"><span class="doc">dump local</span></a>
commands, this body style produces one datum for each of the N
sub-particles in a body particle. The datum has 3 values:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="mi">1</span> <span class="o">=</span> <span class="n">x</span> <span class="n">position</span> <span class="n">of</span> <span class="n">sub</span><span class="o">-</span><span class="n">particle</span>
<span class="mi">2</span> <span class="o">=</span> <span class="n">y</span> <span class="n">position</span> <span class="n">of</span> <span class="n">sub</span><span class="o">-</span><span class="n">particle</span>
<span class="mi">3</span> <span class="o">=</span> <span class="n">z</span> <span class="n">position</span> <span class="n">of</span> <span class="n">sub</span><span class="o">-</span><span class="n">particle</span>
</pre></div>
</div>
<p>These values are the current position of the sub-particle within the
simulation domain, not a displacement from the center-of-mass (COM) of
the body particle itself. These values are calculated using the
current COM and orientation of the body particle.</p>
<p>For images created by the <a class="reference internal" href="dump_image.html"><span class="doc">dump image</span></a> command, if the
<em>body</em> keyword is set, then each body particle is drawn as a
collection of spheres, one for each sub-particle. The size of each
sphere is determined by the <em>bflag1</em> parameter for the <em>body</em> keyword.
The <em>bflag2</em> argument is ignored.</p>
<hr class="docutils" />
<p><strong>Specifics of body style rounded/polygon:</strong></p>
<p>The <em>rounded/polygon</em> body style represents body particles as a convex
polygon with a variable number N &gt; 2 of vertices, which can only be
used for 2d models. One example use of this body style is for 2d
discrete element models, as described in <a class="reference internal" href="#fraige"><span class="std std-ref">Fraige</span></a>. Similar to
body style <em>nparticle</em>, the atom_style body command for this body
style takes two additional arguments:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">atom_style</span> <span class="n">body</span> <span class="n">rounded</span><span class="o">/</span><span class="n">polygon</span> <span class="n">Nmin</span> <span class="n">Nmax</span>
<span class="n">Nmin</span> <span class="o">=</span> <span class="n">minimum</span> <span class="c1"># of vertices in any body in the system</span>
<span class="n">Nmax</span> <span class="o">=</span> <span class="n">maximum</span> <span class="c1"># of vertices in any body in the system</span>
</pre></div>
</div>
<p>The Nmin and Nmax arguments are used to bound the size of data
structures used internally by each particle.</p>
<p>When the <a class="reference internal" href="read_data.html"><span class="doc">read_data</span></a> command reads a data file for this
body style, the following information must be provided for each entry
in the <em>Bodies</em> section of the data file:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">atom</span><span class="o">-</span><span class="n">ID</span> <span class="mi">1</span> <span class="n">M</span>
<span class="n">N</span>
<span class="n">ixx</span> <span class="n">iyy</span> <span class="n">izz</span> <span class="n">ixy</span> <span class="n">ixz</span> <span class="n">iyz</span>
<span class="n">x1</span> <span class="n">y1</span> <span class="n">z1</span>
<span class="o">...</span>
<span class="n">xN</span> <span class="n">yN</span> <span class="n">zN</span>
<span class="n">i</span> <span class="n">j</span> <span class="n">j</span> <span class="n">k</span> <span class="n">k</span> <span class="o">...</span>
<span class="n">radius</span>
</pre></div>
</div>
<p>N is the number of vertices in the body particle. M = 6 + 3*N + 2*N +
1. The integer line has a single value N. The floating point line(s)
list 6 moments of inertia followed by the coordinates of the N
vertices (x1 to zN) as 3N values, followed by 2N vertex indices
corresponding to the end points of the N edges, followed by a single
radius value = the smallest circle encompassing the polygon. That
last value is used to facilitate the body/body contact detection.
These floating-point values can be listed on as many lines as you
wish; see the <a class="reference internal" href="read_data.html"><span class="doc">read_data</span></a> command for more details.</p>
<p>The 6 moments of inertia (ixx,iyy,izz,ixy,ixz,iyz) should be the
values consistent with the current orientation of the rigid body
around its center of mass. The values are with respect to the
simulation box XYZ axes, not with respect to the prinicpal axes of the
rigid body itself. LAMMPS performs the latter calculation internally.
The coordinates of each vertex are specified as its x,y,z displacement
from the center-of-mass of the body particle. The center-of-mass
position of the particle is specified by the x,y,z values in the
<em>Atoms</em> section of the data file.</p>
<p>For example, the following information would specify a square
particles whose edge length is sqrt(2):</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="mi">3</span> <span class="mi">1</span> <span class="mi">27</span>
<span class="mi">4</span>
<span class="mi">1</span> <span class="mi">1</span> <span class="mi">4</span> <span class="mi">0</span> <span class="mi">0</span> <span class="mi">0</span>
<span class="o">-</span><span class="mf">0.7071</span> <span class="o">-</span><span class="mf">0.7071</span> <span class="mi">0</span>
<span class="o">-</span><span class="mf">0.7071</span> <span class="mf">0.7071</span> <span class="mi">0</span>
<span class="mf">0.7071</span> <span class="mf">0.7071</span> <span class="mi">0</span>
<span class="mf">0.7071</span> <span class="o">-</span><span class="mf">0.7071</span> <span class="mi">0</span>
<span class="mi">0</span> <span class="mi">1</span> <span class="mi">1</span> <span class="mi">2</span> <span class="mi">2</span> <span class="mi">3</span> <span class="mi">3</span> <span class="mi">0</span>
<span class="mf">1.0</span>
</pre></div>
</div>
<p>The <span class="xref doc">pair_style body/rounded/polygon</span> command
can be used with this body style to compute body/body interactions.</p>
<p>For output purposes via the <a class="reference internal" href="compute_body_local.html"><span class="doc">compute body/local</span></a> and <a class="reference internal" href="dump.html"><span class="doc">dump local</span></a>
commands, this body style produces one datum for each of the N
sub-particles in a body particle. The datum has 3 values:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="mi">1</span> <span class="o">=</span> <span class="n">x</span> <span class="n">position</span> <span class="n">of</span> <span class="n">vertex</span>
<span class="mi">2</span> <span class="o">=</span> <span class="n">y</span> <span class="n">position</span> <span class="n">of</span> <span class="n">vertex</span>
<span class="mi">3</span> <span class="o">=</span> <span class="n">z</span> <span class="n">position</span> <span class="n">of</span> <span class="n">vertex</span>
</pre></div>
</div>
<p>These values are the current position of the vertex within the
simulation domain, not a displacement from the center-of-mass (COM) of
the body particle itself. These values are calculated using the
current COM and orientation of the body particle.</p>
<p>For images created by the <a class="reference internal" href="dump_image.html"><span class="doc">dump image</span></a> command, if the
<em>body</em> keyword is set, then each body particle is drawn as a convex
polygon consisting of N line segments. Note that the line segments
are drawn between the N vertices, which does not correspond exactly to
the physical extent of the body (because the <a class="reference external" href="pair_body_rounded_polygon.cpp">pair_style rounded/polygon</a> defines finite-size
spheres at those point and the line segments between the spheres are
tangent to the spheres). The drawn diameter of each line segment is
determined by the <em>bflag1</em> parameter for the <em>body</em> keyword. The
<em>bflag2</em> argument is ignored.</p>
<hr class="docutils" />
<p id="fraige"><strong>(Fraige)</strong> F. Y. Fraige, P. A. Langston, A. J. Matchett, J. Dodds,
Particuology, 6, 455 (2008).</p>
</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>