Updated compute hexorder/atom, added compute orientorder/atom

git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@14251 f3b2605a-c512-4ea7-a41b-209d697bcdaa
This commit is contained in:
athomps 2015-11-15 22:34:00 +00:00
parent da1a3ac83a
commit 70aba85d31
11 changed files with 936 additions and 335 deletions

Binary file not shown.

Before

Width:  |  Height:  |  Size: 6.4 KiB

After

Width:  |  Height:  |  Size: 48 KiB

View File

@ -2,7 +2,7 @@
\begin{document}
$$
q_n = \frac{1}{n}\sum_{j = 1}^{n} e^{n i \theta({\bf r}_{ij})}
q_n = \frac{1}{nnn}\sum_{j = 1}^{nnn} e^{n i \theta({\bf r}_{ij})}
$$
\end{document}

BIN
doc/Eqs/orientorder.jpg Normal file

Binary file not shown.

After

Width:  |  Height:  |  Size: 125 KiB

12
doc/Eqs/orientorder.tex Normal file
View File

@ -0,0 +1,12 @@
\documentclass[12pt]{article}
\begin{document}
$$
\bar{Y}_{lm} = \frac{1}{nnn}\sum_{j = 1}^{nnn} Y_{lm}( \theta( {\bf r}_{ij} ), \phi( {\bf r}_{ij} ) )
$$
$$
Q^2_l = \frac{4 \pi}{2 l + 1} \sum_{m = -l}^{m = l} \bar{Y}_{lm} \bar{Y}^*_{lm}
$$
\end{document}

View File

@ -1,286 +1,125 @@
<HTML>
<CENTER><A HREF = "http://lammps.sandia.gov">LAMMPS WWW Site</A> - <A HREF = "Manual.html">LAMMPS Documentation</A> - <A HREF = "Section_commands.html#comm">LAMMPS Commands</A>
</CENTER>
<!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>compute hexorder/atom command &mdash; LAMMPS 15 May 2015 version 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 15 May 2015 version 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>compute hexorder/atom 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="compute-hexorder-atom-command">
<span id="index-0"></span><h1>compute hexorder/atom command<a class="headerlink" href="#compute-hexorder-atom-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>compute ID group-ID hexorder/atom keyword values ...
</pre></div>
</div>
<ul class="simple">
<li>ID, group-ID are documented in <a class="reference internal" href="compute.html"><em>compute</em></a> command</li>
<li>hexorder/atom = style name of this compute command</li>
<li>zero or more keyword/value pairs may be appended</li>
<li>keyword = <em>degree</em></li>
</ul>
<pre class="literal-block">
<em>n</em> value = degree of order parameter
</pre>
</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>compute 1 all hexorder/atom
compute 1 all hexorder/atom n 4
</pre></div>
</div>
</div>
<div class="section" id="description">
<h2>Description<a class="headerlink" href="#description" title="Permalink to this headline"></a></h2>
<p>Define a computation that calculates <em>qn</em> the bond-orientational
order parameter for each atom in a group. The hexatic (<em>n</em> = 6) order
parameter was introduced by <a class="reference internal" href="#nelson"><span>Nelson and Halperin</span></a> as a way to detect
hexagonal symmetry in two-dimensional systems. For each atom, <em>qn</em>
is a complex number (stored as two real numbers) defined as follows:</p>
<img alt="_images/hexorder.jpg" class="align-center" src="_images/hexorder.jpg" />
<p>where the sum is over the <em>n</em> nearest neighbors
<HR>
<H3>compute hexorder/atom command
</H3>
<P><B>Syntax:</B>
</P>
<PRE>compute ID group-ID hexorder/atom keyword values ...
</PRE>
<UL><LI>ID, group-ID are documented in <A HREF = "compute.html">compute</A> command
<LI>hexorder/atom = style name of this compute command
<LI>one or more keyword/value pairs may be appended
<PRE>keyword = <I>n</I> or <I>nnn</I> or <I>cutoff</I>
<I>cutoff</I> value = distance cutoff
<I>nnn</I> value = number of nearest neighbors
<I>n</I> value = degree of order parameter
</PRE>
</UL>
<P><B>Examples:</B>
</P>
<PRE>compute 1 all hexorder/atom
compute 1 all hexorder/atom n 4 nnn 4 cutoff 1.2
</PRE>
<P><B>Description:</B>
</P>
<P>Define a computation that calculates <I>qn</I> the bond-orientational
order parameter for each atom in a group. The hexatic (<I>n</I> = 6) order
parameter was introduced by <A HREF = "#Nelson">Nelson and Halperin</A> as a way to detect
hexagonal symmetry in two-dimensional systems. For each atom, <I>qn</I>
is a complex number (stored as two real numbers) defined as follows:
</P>
<CENTER><IMG SRC = "Eqs/hexorder.jpg">
</CENTER>
<P>where the sum is over the <I>nnn</I> nearest neighbors
of the central atom. The angle theta
is formed by the bond vector rij and the <em>x</em> axis. theta is calculated
only using the <em>x</em> and <em>y</em> components, whereas the distance from the
central atom is calculated using all three
<em>x</em>, <em>y</em>, and <em>z</em> components of the bond vector.
Neighbor atoms not in the group
are included in the order parameter of atoms in the group.</p>
<p>The optional keyword <em>n</em> sets the degree of the order parameter.
The default value is 6. For a perfect hexagonal lattice,
<a href="#id1"><span class="problematic" id="id2">*</span></a>q*6 = exp(6 i phi) for all atoms, where the constant 0 &lt; phi &lt; pi/3
depends only on the orientation of the lattice relative to the x axis.
In an isotropic liquid, local neighborhoods may still exhibit
is formed by the bond vector rij and the <I>x</I> axis. theta is calculated
only using the <I>x</I> and <I>y</I> components, whereas the distance from the
central atom is calculated using all three
<I>x</I>, <I>y</I>, and <I>z</I> components of the bond vector.
Neighbor atoms not in the group
are included in the order parameter of atoms in the group.
</P>
<P>The optional keyword <I>cutoff</I> defines the distance cutoff
used when searching for neighbors. The default value, also
the maximum allowable value, is the cutoff specified
by the pair style.
</P>
<P>The optional keyword <I>nnn</I> defines the number of nearest
neighbors used to calculate <I>qn</I>. The default value is 6.
If the value is NULL, then all neighbors up to the
distance cutoff are used.
</P>
<P>The optional keyword <I>n</I> sets the degree of the order parameter.
The default value is 6. For a perfect hexagonal lattice with
<I>nnn</I> = 6,
<I>q</I>6 = exp(6 i phi) for all atoms, where the constant 0 < phi < pi/3
depends only on the orientation of the lattice relative to the <I>x</I> axis.
In an isotropic liquid, local neighborhoods may still exhibit
weak hexagonal symmetry, but because the orientational correlation
decays quickly with distance, the value of phi will be different for
different atoms, and <a href="#id3"><span class="problematic" id="id4">|&lt;*q*6&gt;|</span></a> &lt;&lt; 1.</p>
<p>The value of <em>qn</em> will be zero for atoms not in the
specified compute group. If the atom has less than <em>n</em> neighbors (within
the potential cutoff), then <em>qn</em> is set to zero.</p>
<p>The neighbor list needed to compute this quantity is constructed each
different atoms, and so when <I>q</I>6 is averaged over all the atoms
in the system, |<<I>q</I>6>| << 1.
</P>
<P>The value of <I>qn</I> is set to zero for atoms not in the
specified compute group, as well as for atoms that have less than
<I>nnn</I> neighbors within the distance cutoff.
</P>
<P>The neighbor list needed to compute this quantity is constructed each
time the calculation is performed (i.e. each time a snapshot of atoms
is dumped). Thus it can be inefficient to compute/dump this quantity
too frequently.</p>
<div class="admonition warning">
<p class="first admonition-title">Warning</p>
<p class="last">If you have a bonded system, then the settings of
<a class="reference internal" href="special_bonds.html"><em>special_bonds</em></a> command can remove pairwise
too frequently.
</P>
<P>IMPORTANT NOTE: If you have a bonded system, then the settings of
<A HREF = "special_bonds.html">special_bonds</A> command can remove pairwise
interactions between atoms in the same bond, angle, or dihedral. This
is the default setting for the <a class="reference internal" href="special_bonds.html"><em>special_bonds</em></a>
is the default setting for the <A HREF = "special_bonds.html">special_bonds</A>
command, and means those pairwise interactions do not appear in the
neighbor list. Because this fix uses the neighbor list, it also means
those pairs will not be included in the order parameter. One way
to get around this, is to write a dump file, and use the
<a class="reference internal" href="rerun.html"><em>rerun</em></a> command to compute the order parameter for snapshots
those pairs will not be included in the order parameter. This difficulty
can be circumvented by writing a dump file, and using the
<A HREF = "rerun.html">rerun</A> command to compute the order parameter for snapshots
in the dump file. The rerun script can use a
<a class="reference internal" href="special_bonds.html"><em>special_bonds</em></a> command that includes all pairs in
the neighbor list.</p>
</div>
<p><strong>Output info:</strong></p>
<p>This compute calculates a per-atom array with 2 columns, giving the
real and imaginary parts of <em>qn</em>, respectively.</p>
<p>These values can be accessed by any command that uses
per-atom values from a compute as input. See <a class="reference internal" href="Section_howto.html#howto-15"><span>Section_howto 15</span></a> for an overview of LAMMPS output
options.</p>
<p>The per-atom array contain pairs of numbers representing the
real and imaginary parts of <em>qn</em>, a complex number subject to the
constraint <a href="#id5"><span class="problematic" id="id6">|*qn*|</span></a> &lt;= 1.</p>
</div>
<div class="section" id="restrictions">
<h2>Restrictions<a class="headerlink" href="#restrictions" title="Permalink to this headline"></a></h2>
<blockquote>
<div>none</div></blockquote>
</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="compute_coord_atom.html"><em>compute coord/atom</em></a>, <a class="reference internal" href="compute_centro_atom.html"><em>compute centro/atom</em></a></p>
</div>
<div class="section" id="default">
<h2>Default<a class="headerlink" href="#default" title="Permalink to this headline"></a></h2>
<p>The option default is <em>n</em> = 6.</p>
<hr class="docutils" />
<p id="nelson"><strong>(Nelson)</strong> Nelson, Halperin, Phys Rev B, 19, 2457 (1979).</p>
</div>
</div>
<A HREF = "special_bonds.html">special_bonds</A> command that includes all pairs in
the neighbor list.
</P>
<P><B>Output info:</B>
</P>
<P>This compute calculates a per-atom array with 2 columns, giving the
real and imaginary parts <I>qn</I>, a complex number restricted to the
unit disk of the complex plane i.e. Re(<I>qn</I>)^2 + Im(<I>qn</I>)^2 <= 1 .
</P>
<P>These values can be accessed by any command that uses
per-atom values from a compute as input. See <A HREF = "Section_howto.html#howto_15">Section_howto
15</A> for an overview of LAMMPS output
options.
</P>
<P><B>Restrictions:</B> none
</P>
<P><B>Related commands:</B>
</P>
<P><A HREF = "compute_orientorder_atom.html">compute orientorder/atom</A>, <A HREF = "compute_coord_atom.html">compute coord/atom</A>, <A HREF = "compute_centro_atom.html">compute centro/atom</A>
</P>
<P><B>Default:</B>
</P>
<P>The option defaults are <I>n</I> = 6, <I>nnn</I> = 6, <I>cutoff</I> = pair style cutoff
</P>
<HR>
<A NAME = "Nelson"></A>
</div>
</div>
<footer>
<hr/>
<div role="contentinfo">
<p>
&copy; Copyright .
</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:'15 May 2015 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>
<P><B>(Nelson)</B> Nelson, Halperin, Phys Rev B, 19, 2457 (1979).
</P>
</HTML>

View File

@ -14,15 +14,18 @@ compute ID group-ID hexorder/atom keyword values ... :pre
ID, group-ID are documented in "compute"_compute.html command :ulb,l
hexorder/atom = style name of this compute command :l
zero or more keyword/value pairs may be appended :l
keyword = {degree} :l
one or more keyword/value pairs may be appended :l
keyword = {n} or {nnn} or {cutoff}
{cutoff} value = distance cutoff
{nnn} value = number of nearest neighbors
{n} value = degree of order parameter :pre
:ule
[Examples:]
compute 1 all hexorder/atom
compute 1 all hexorder/atom n 4 :pre
compute 1 all hexorder/atom n 4 nnn 4 cutoff 1.2 :pre
[Description:]
@ -34,7 +37,7 @@ is a complex number (stored as two real numbers) defined as follows:
:c,image(Eqs/hexorder.jpg)
where the sum is over the {n} nearest neighbors
where the sum is over the {nnn} nearest neighbors
of the central atom. The angle theta
is formed by the bond vector rij and the {x} axis. theta is calculated
only using the {x} and {y} components, whereas the distance from the
@ -43,18 +46,30 @@ central atom is calculated using all three
Neighbor atoms not in the group
are included in the order parameter of atoms in the group.
The optional keyword {cutoff} defines the distance cutoff
used when searching for neighbors. The default value, also
the maximum allowable value, is the cutoff specified
by the pair style.
The optional keyword {nnn} defines the number of nearest
neighbors used to calculate {qn}. The default value is 6.
If the value is NULL, then all neighbors up to the
distance cutoff are used.
The optional keyword {n} sets the degree of the order parameter.
The default value is 6. For a perfect hexagonal lattice,
The default value is 6. For a perfect hexagonal lattice with
{nnn} = 6,
{q}6 = exp(6 i phi) for all atoms, where the constant 0 < phi < pi/3
depends only on the orientation of the lattice relative to the x axis.
depends only on the orientation of the lattice relative to the {x} axis.
In an isotropic liquid, local neighborhoods may still exhibit
weak hexagonal symmetry, but because the orientational correlation
decays quickly with distance, the value of phi will be different for
different atoms, and |<{q}6>| << 1.
different atoms, and so when {q}6 is averaged over all the atoms
in the system, |<{q}6>| << 1.
The value of {qn} will be zero for atoms not in the
specified compute group. If the atom has less than {n} neighbors (within
the potential cutoff), then {qn} is set to zero.
The value of {qn} is set to zero for atoms not in the
specified compute group, as well as for atoms that have less than
{nnn} neighbors within the distance cutoff.
The neighbor list needed to compute this quantity is constructed each
time the calculation is performed (i.e. each time a snapshot of atoms
@ -67,8 +82,8 @@ interactions between atoms in the same bond, angle, or dihedral. This
is the default setting for the "special_bonds"_special_bonds.html
command, and means those pairwise interactions do not appear in the
neighbor list. Because this fix uses the neighbor list, it also means
those pairs will not be included in the order parameter. One way
to get around this, is to write a dump file, and use the
those pairs will not be included in the order parameter. This difficulty
can be circumvented by writing a dump file, and using the
"rerun"_rerun.html command to compute the order parameter for snapshots
in the dump file. The rerun script can use a
"special_bonds"_special_bonds.html command that includes all pairs in
@ -77,26 +92,23 @@ the neighbor list.
[Output info:]
This compute calculates a per-atom array with 2 columns, giving the
real and imaginary parts of {qn}, respectively.
real and imaginary parts {qn}, a complex number restricted to the
unit disk of the complex plane i.e. Re({qn})^2 + Im({qn})^2 <= 1 .
These values can be accessed by any command that uses
per-atom values from a compute as input. See "Section_howto
15"_Section_howto.html#howto_15 for an overview of LAMMPS output
options.
The per-atom array contain pairs of numbers representing the
real and imaginary parts of {qn}, a complex number subject to the
constraint |{qn}| <= 1.
[Restrictions:] none
[Related commands:]
"compute coord/atom"_compute_coord_atom.html, "compute centro/atom"_compute_centro_atom.html
"compute orientorder/atom"_compute_orientorder_atom.html, "compute coord/atom"_compute_coord_atom.html, "compute centro/atom"_compute_centro_atom.html
[Default:]
The option default is {n} = 6.
The option defaults are {n} = 6, {nnn} = 6, {cutoff} = pair style cutoff
:line

View File

@ -0,0 +1,115 @@
"LAMMPS WWW Site"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c
:link(lws,http://lammps.sandia.gov)
:link(ld,Manual.html)
:link(lc,Section_commands.html#comm)
:line
compute orientorder/atom command :h3
[Syntax:]
compute ID group-ID orientorder/atom keyword values ... :pre
ID, group-ID are documented in "compute"_compute.html command :ulb,l
orientorder/atom = style name of this compute command :l
one or more keyword/value pairs may be appended :l
keyword = {cutoff} or {nnn} or {ql}
{cutoff} value = distance cutoff
{nnn} value = number of nearest neighbors
{degrees} values = nlvalues, l1, l2,... :pre
:ule
[Examples:]
compute 1 all orientorder/atom
compute 1 all orientorder/atom degrees 5 4 6 8 10 12 nnn NULL cutoff 1.5 :pre
[Description:]
Define a computation that calculates a set of bond-orientational
order parameters {Ql} for each atom in a group. These order parameters
were introduced by "Steinhardt et al."_#Steinhardt as a way to
characterize the local orientational order in atomic structures.
For each atom, {Ql} is a real number defined as follows:
:c,image(Eqs/orientorder.jpg)
The first equation defines the spherical harmonic order parameters.
These are complex number components of the 3D analog of the 2D order
parameter {qn}, which is implemented as LAMMPS compute
"hexorder/atom"_compute_hexorder_atom.html.
The summation is over the {nnn} nearest
neighbors of the central atom.
The angles theta and phi are the standard spherical polar angles
defining the direction of the bond vector {rij}.
The second equation defines the square power of {Ql}, which are
rotationally invariant scalar quantities obtained by summing
over all the components of degree {l}.
The optional keyword {cutoff} defines the distance cutoff
used when searching for neighbors. The default value, also
the maximum allowable value, is the cutoff specified
by the pair style.
The optional keyword {nnn} defines the number of nearest
neighbors used to calculate {Ql}. The default value is 12.
If the value is NULL, then all neighbors up to the
specified distance cutoff are used.
The optional keyword {degrees} defines the list of order parameters to
be computed. The first argument {nlvalues} is the number of order
parameters. This is followed by that number of integers giving the
degree of each order parameter. Because {Q}2 and all odd-degree
order parameters are zero for atoms in cubic crystals
(see "Steinhardt"_#Steinhardt), the default order parameters
are {Q}4, {Q}6, {Q}8, {Q}10, and {Q}12.
The value of {Ql} is set to zero for atoms not in the
specified compute group, as well as for atoms that have less than
{nnn} neighbors within the distance cutoff.
The neighbor list needed to compute this quantity is constructed each
time the calculation is performed (i.e. each time a snapshot of atoms
is dumped). Thus it can be inefficient to compute/dump this quantity
too frequently.
IMPORTANT NOTE: If you have a bonded system, then the settings of
"special_bonds"_special_bonds.html command can remove pairwise
interactions between atoms in the same bond, angle, or dihedral. This
is the default setting for the "special_bonds"_special_bonds.html
command, and means those pairwise interactions do not appear in the
neighbor list. Because this fix uses the neighbor list, it also means
those pairs will not be included in the order parameter. This difficulty
can be circumvented by writing a dump file, and using the
"rerun"_rerun.html command to compute the order parameter for snapshots
in the dump file. The rerun script can use a
"special_bonds"_special_bonds.html command that includes all pairs in
the neighbor list.
[Output info:]
This compute calculates a per-atom array with {nlvalues} columns, giving the
{Ql} values for each atom, which are real numbers on the range 0 <= {Ql} <= 1.
These values can be accessed by any command that uses
per-atom values from a compute as input. See "Section_howto
15"_Section_howto.html#howto_15 for an overview of LAMMPS output
options.
[Restrictions:] none
[Related commands:]
"compute coord/atom"_compute_coord_atom.html, "compute centro/atom"_compute_centro_atom.html, "compute hexorder/atom"_compute_hexorder_atom.html
[Default:]
The option defaults are {cutoff} = pair style cutoff, {nnn} = 12, {degrees} = 5 4 6 8 9 10 12 i.e. {Q}4, {Q}6, {Q}8, {Q}10, and {Q}12.
:line
:link(Steinhardt)
[(Steinhardt)] P. Steinhardt, D. Nelson, and M. Ronchetti, Phys. Rev. B 28, 784 (1983).

View File

@ -15,7 +15,6 @@
Contributing author: Aidan Thompson (SNL)
------------------------------------------------------------------------- */
#include <math.h>
#include <complex>
#include <string.h>
#include <stdlib.h>
@ -31,8 +30,16 @@
#include "comm.h"
#include "memory.h"
#include "error.h"
#include "math_const.h"
#ifdef DBL_EPSILON
#define MY_EPSILON (10.0*DBL_EPSILON)
#else
#define MY_EPSILON (10.0*2.220446049250313e-16)
#endif
using namespace LAMMPS_NS;
using namespace MathConst;
/* ---------------------------------------------------------------------- */
@ -41,18 +48,37 @@ ComputeHexOrderAtom::ComputeHexOrderAtom(LAMMPS *lmp, int narg, char **arg) :
{
if (narg < 3 ) error->all(FLERR,"Illegal compute hexorder/atom command");
ndegree = 6;
nnn = 6;
cutsq = 0.0;
// process optional args
int iarg = 3;
while (iarg < narg) {
if (strcmp(arg[iarg],"n") == 0) {
if (iarg+1 > narg) error->all(FLERR,"Illegal compute hexorder/atom command");
nnn = force->numeric(FLERR,arg[iarg+1]);
if (nnn < 0)
if (iarg+2 > narg) error->all(FLERR,"Illegal compute hexorder/atom command");
ndegree = force->numeric(FLERR,arg[iarg+1]);
if (ndegree < 0)
error->all(FLERR,"Illegal compute hexorder/atom command");
iarg += 2;
} else if (strcmp(arg[iarg],"nnn") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal compute hexorder/atom command");
if (strcmp(arg[iarg+1],"NULL") == 0)
nnn = 0;
else {
nnn = force->numeric(FLERR,arg[iarg+1]);
if (nnn < 0)
error->all(FLERR,"Illegal compute hexorder/atom command");
}
iarg += 2;
} else if (strcmp(arg[iarg],"cutoff") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal compute hexorder/atom command");
double cutoff = force->numeric(FLERR,arg[iarg+1]);
if (cutoff <= 0.0)
error->all(FLERR,"Illegal compute hexorder/atom command");
cutsq = cutoff*cutoff;
iarg += 2;
} else error->all(FLERR,"Illegal compute hexorder/atom command");
}
@ -61,7 +87,7 @@ ComputeHexOrderAtom::ComputeHexOrderAtom(LAMMPS *lmp, int narg, char **arg) :
size_peratom_cols = ncol;
nmax = 0;
q6array = NULL;
qnarray = NULL;
maxneigh = 0;
distsq = NULL;
nearest = NULL;
@ -71,7 +97,7 @@ ComputeHexOrderAtom::ComputeHexOrderAtom(LAMMPS *lmp, int narg, char **arg) :
ComputeHexOrderAtom::~ComputeHexOrderAtom()
{
memory->destroy(q6array);
memory->destroy(qnarray);
memory->destroy(distsq);
memory->destroy(nearest);
}
@ -82,6 +108,10 @@ void ComputeHexOrderAtom::init()
{
if (force->pair == NULL)
error->all(FLERR,"Compute hexorder/atom requires a pair style be defined");
if (cutsq == 0.0) cutsq = force->pair->cutforce * force->pair->cutforce;
else if (sqrt(cutsq) > force->pair->cutforce)
error->all(FLERR,
"Compute hexorder/atom cutoff is longer than pairwise cutoff");
// need an occasional full neighbor list
@ -119,10 +149,10 @@ void ComputeHexOrderAtom::compute_peratom()
// grow order parameter array if necessary
if (atom->nlocal > nmax) {
memory->destroy(q6array);
memory->destroy(qnarray);
nmax = atom->nmax;
memory->create(q6array,nmax,ncol,"hexorder/atom:q6array");
array_atom = q6array;
memory->create(qnarray,nmax,ncol,"hexorder/atom:qnarray");
array_atom = qnarray;
}
// invoke full neighbor list (will copy or build if necessary)
@ -139,11 +169,10 @@ void ComputeHexOrderAtom::compute_peratom()
double **x = atom->x;
int *mask = atom->mask;
double cutsq = force->pair->cutforce * force->pair->cutforce;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
double* q6 = q6array[i];
double* qn = qnarray[i];
if (mask[i] & groupbit) {
xtmp = x[i][0];
ytmp = x[i][1];
@ -157,8 +186,8 @@ void ComputeHexOrderAtom::compute_peratom()
memory->destroy(distsq);
memory->destroy(nearest);
maxneigh = jnum;
memory->create(distsq,maxneigh,"hexcoord/atom:distsq");
memory->create(nearest,maxneigh,"hexcoord/atom:nearest");
memory->create(distsq,maxneigh,"hexorder/atom:distsq");
memory->create(nearest,maxneigh,"hexorder/atom:nearest");
}
// loop over list of all neighbors within force cutoff
@ -183,62 +212,63 @@ void ComputeHexOrderAtom::compute_peratom()
// if not nnn neighbors, order parameter = 0;
if (ncount < nnn) {
q6[0] = q6[1] = 0.0;
qn[0] = qn[1] = 0.0;
continue;
}
// store nnn nearest neighs in 1st nnn locations of distsq and nearest
// if nnn > 0, use only nearest nnn neighbors
select2(nnn,ncount,distsq,nearest);
if (nnn > 0) {
select2(nnn,ncount,distsq,nearest);
ncount = nnn;
}
double usum = 0.0;
double vsum = 0.0;
for (jj = 0; jj < nnn; jj++) {
for (jj = 0; jj < ncount; jj++) {
j = nearest[jj];
j &= NEIGHMASK;
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
double u, v;
calc_qn(delx, dely, u, v);
calc_qn_complex(delx, dely, u, v);
usum += u;
vsum += v;
}
q6[0] = usum/nnn;
q6[1] = vsum/nnn;
qn[0] = usum/nnn;
qn[1] = vsum/nnn;
}
}
}
// this might be faster than pow(std::complex) on some platforms
// calculate order parameter using std::complex::pow function
inline void ComputeHexOrderAtom::calc_q6(double delx, double dely, double &u, double &v) {
inline void ComputeHexOrderAtom::calc_qn_complex(double delx, double dely, double &u, double &v) {
double rinv = 1.0/sqrt(delx*delx+dely*dely);
double x = delx*rinv;
double y = dely*rinv;
double a = x*x;
double b1 = y*y;
double b2 = b1*b1;
double b3 = b2*b1;
// (x + i y)^6 coeffs: 1, 6, -15, -20, 15, 6, -1
u = (( a - 15*b1)*a + 15*b2)*a - b3;
v = ((6*a - 20*b1)*a + 6*b2)*x*y;
}
inline void ComputeHexOrderAtom::calc_qn(double delx, double dely, double &u, double &v) {
double rinv = 1.0/sqrt(delx*delx+dely*dely);
double x = delx*rinv;
double y = dely*rinv;
// std::complex<double> z = std::complex<double>(x, y);
std::complex<double> z(x, y);
std::complex<double> zn = pow(z,nnn);
std::complex<double> zn = pow(z, nnn);
u = real(zn);
v = imag(zn);
}
// calculate order parameter using trig functions
// this is usually slower, but can be used if <complex> not available
inline void ComputeHexOrderAtom::calc_qn_trig(double delx, double dely, double &u, double &v) {
double ntheta;
if(fabs(delx) <= MY_EPSILON) {
if(dely > 0.0) ntheta = nnn * MY_PI / 2.0;
else ntheta = nnn * 3.0 * MY_PI / 2.0;
} else ntheta = nnn * atan(dely / delx);
ntheta = nnn * atan(dely / delx);
u = cos(ntheta);
v = sin(ntheta);
}
/* ----------------------------------------------------------------------
select2 routine from Numerical Recipes (slightly modified)
find k smallest values in array of length n
@ -310,5 +340,8 @@ void ComputeHexOrderAtom::select2(int k, int n, double *arr, int *iarr)
double ComputeHexOrderAtom::memory_usage()
{
double bytes = ncol*nmax * sizeof(double);
bytes += maxneigh * sizeof(double);
bytes += maxneigh * sizeof(int);
return bytes;
}

View File

@ -34,16 +34,16 @@ class ComputeHexOrderAtom : public Compute {
double memory_usage();
private:
int nmax,maxneigh,ncol,nnn;
int nmax,maxneigh,ncol,nnn,ndegree;
double cutsq;
class NeighList *list;
double *distsq;
int *nearest;
double **q6array;
double **qnarray;
void calc_q6(double, double, double&, double&);
void calc_q4(double, double, double&, double&);
void calc_qn(double, double, double&, double&);
void calc_qn_complex(double, double, double&, double&);
void calc_qn_trig(double, double, double&, double&);
void select2(int, int, double *, int *);
};

View File

@ -0,0 +1,510 @@
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
/* ----------------------------------------------------------------------
Contributing author: Aidan Thompson (SNL)
Axel Kohlmeyer (Temple U)
------------------------------------------------------------------------- */
#include <string.h>
#include <stdlib.h>
#include "compute_orientorder_atom.h"
#include "atom.h"
#include "update.h"
#include "modify.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "neigh_request.h"
#include "force.h"
#include "pair.h"
#include "comm.h"
#include "memory.h"
#include "error.h"
#include "math_const.h"
using namespace LAMMPS_NS;
using namespace MathConst;
using namespace std;
#ifdef DBL_EPSILON
#define MY_EPSILON (10.0*DBL_EPSILON)
#else
#define MY_EPSILON (10.0*2.220446049250313e-16)
#endif
/* ---------------------------------------------------------------------- */
ComputeOrientOrderAtom::ComputeOrientOrderAtom(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg < 3 ) error->all(FLERR,"Illegal compute orientorder/atom command");
// set default values for optional args
nnn = 12;
cutsq = 0.0;
// specify which orders to request
nqlist = 5;
memory->create(qlist,nqlist,"orientorder/atom:qlist");
qlist[0] = 4;
qlist[1] = 6;
qlist[2] = 8;
qlist[3] = 10;
qlist[4] = 12;
qmax = 12;
// process optional args
int iarg = 3;
while (iarg < narg) {
if (strcmp(arg[iarg],"nnn") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
if (strcmp(arg[iarg+1],"NULL") == 0)
nnn = 0;
else {
nnn = force->numeric(FLERR,arg[iarg+1]);
if (nnn <= 0)
error->all(FLERR,"Illegal compute orientorder/atom command");
}
iarg += 2;
} else if (strcmp(arg[iarg],"degrees") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
nqlist = force->numeric(FLERR,arg[iarg+1]);
if (nqlist <= 0) error->all(FLERR,"Illegal compute orientorder/atom command");
memory->destroy(qlist);
memory->create(qlist,nqlist,"orientorder/atom:qlist");
iarg += 2;
if (iarg+nqlist > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
qmax = 0;
for (int iw = 0; iw < nqlist; iw++) {
qlist[iw] = force->numeric(FLERR,arg[iarg+iw]);
if (qlist[iw] < 0)
error->all(FLERR,"Illegal compute orientorder/atom command");
if (qlist[iw] > qmax) qmax = qlist[iw];
}
iarg += nqlist;
} else if (strcmp(arg[iarg],"cutoff") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
double cutoff = force->numeric(FLERR,arg[iarg+1]);
if (cutoff <= 0.0) error->all(FLERR,"Illegal compute orientorder/atom command");
cutsq = cutoff*cutoff;
iarg += 2;
} else error->all(FLERR,"Illegal compute orientorder/atom command");
}
ncol = nqlist;
peratom_flag = 1;
size_peratom_cols = ncol;
nmax = 0;
qnarray = NULL;
maxneigh = 0;
distsq = NULL;
nearest = NULL;
rlist = NULL;
qnm_r = NULL;
qnm_i = NULL;
}
/* ---------------------------------------------------------------------- */
ComputeOrientOrderAtom::~ComputeOrientOrderAtom()
{
memory->destroy(qnarray);
memory->destroy(distsq);
memory->destroy(nearest);
memory->destroy(qlist);
memory->destroy(qnm_r);
memory->destroy(qnm_i);
}
/* ---------------------------------------------------------------------- */
void ComputeOrientOrderAtom::init()
{
if (force->pair == NULL)
error->all(FLERR,"Compute orientorder/atom requires a pair style be defined");
if (cutsq == 0.0) cutsq = force->pair->cutforce * force->pair->cutforce;
else if (sqrt(cutsq) > force->pair->cutforce)
error->all(FLERR,
"Compute orientorder/atom cutoff is longer than pairwise cutoff");
memory->create(qnm_r,qmax,2*qmax+1,"orientorder/atom:qnm_r");
memory->create(qnm_i,qmax,2*qmax+1,"orientorder/atom:qnm_i");
// need an occasional full neighbor list
int irequest = neighbor->request(this,instance_me);
neighbor->requests[irequest]->pair = 0;
neighbor->requests[irequest]->compute = 1;
neighbor->requests[irequest]->half = 0;
neighbor->requests[irequest]->full = 1;
neighbor->requests[irequest]->occasional = 1;
int count = 0;
for (int i = 0; i < modify->ncompute; i++)
if (strcmp(modify->compute[i]->style,"orientorder/atom") == 0) count++;
if (count > 1 && comm->me == 0)
error->warning(FLERR,"More than one compute orientorder/atom");
}
/* ---------------------------------------------------------------------- */
void ComputeOrientOrderAtom::init_list(int id, NeighList *ptr)
{
list = ptr;
}
/* ---------------------------------------------------------------------- */
void ComputeOrientOrderAtom::compute_peratom()
{
int i,j,m,ii,jj,inum,jnum;
double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
int *ilist,*jlist,*numneigh,**firstneigh;
invoked_peratom = update->ntimestep;
// grow order parameter array if necessary
if (atom->nlocal > nmax) {
memory->destroy(qnarray);
nmax = atom->nmax;
memory->create(qnarray,nmax,ncol,"orientorder/atom:qnarray");
array_atom = qnarray;
}
// invoke full neighbor list (will copy or build if necessary)
neighbor->build_one(list);
inum = list->inum;
ilist = list->ilist;
numneigh = list->numneigh;
firstneigh = list->firstneigh;
// compute order parameter for each atom in group
// use full neighbor list to count atoms less than cutoff
double **x = atom->x;
int *mask = atom->mask;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
double* qn = qnarray[i];
if (mask[i] & groupbit) {
xtmp = x[i][0];
ytmp = x[i][1];
ztmp = x[i][2];
jlist = firstneigh[i];
jnum = numneigh[i];
// insure distsq and nearest arrays are long enough
if (jnum > maxneigh) {
memory->destroy(distsq);
memory->destroy(rlist);
memory->destroy(nearest);
maxneigh = jnum;
memory->create(distsq,maxneigh,"orientorder/atom:distsq");
memory->create(rlist,maxneigh,3,"orientorder/atom:rlist");
memory->create(nearest,maxneigh,"orientorder/atom:nearest");
}
// loop over list of all neighbors within force cutoff
// distsq[] = distance sq to each
// rlist[] = distance vector to each
// nearest[] = atom indices of neighbors
int ncount = 0;
for (jj = 0; jj < jnum; jj++) {
j = jlist[jj];
j &= NEIGHMASK;
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
delz = ztmp - x[j][2];
rsq = delx*delx + dely*dely + delz*delz;
if (rsq < cutsq) {
distsq[ncount] = rsq;
rlist[ncount][0] = delx;
rlist[ncount][1] = dely;
rlist[ncount][2] = delz;
nearest[ncount++] = j;
}
}
// if not nnn neighbors, order parameter = 0;
if (ncount < nnn) {
for (int iw = 0; iw < nqlist; iw++)
qn[iw] = 0.0;
continue;
}
// if nnn > 0, use only nearest nnn neighbors
if (nnn > 0) {
select3(nnn,ncount,distsq,nearest,rlist);
calc_boop(rlist, nnn, qn, qlist, nqlist);
} else
calc_boop(rlist, ncount, qn, qlist, nqlist);
}
}
}
/* ----------------------------------------------------------------------
memory usage of local atom-based array
------------------------------------------------------------------------- */
double ComputeOrientOrderAtom::memory_usage()
{
double bytes = ncol*nmax * sizeof(double);
bytes += (qmax*(2*qmax+1)+maxneigh*4) * sizeof(double);
bytes += (nqlist+maxneigh) * sizeof(int);
return bytes;
}
/* ----------------------------------------------------------------------
select3 routine from Numerical Recipes (slightly modified)
find k smallest values in array of length n
sort auxiliary arrays at same time
------------------------------------------------------------------------- */
// Use no-op do while to create single statement
#define SWAP(a,b) do { \
tmp = a; a = b; b = tmp; \
} while(0)
#define ISWAP(a,b) do { \
itmp = a; a = b; b = itmp; \
} while(0)
#define SWAP3(a,b) do { \
tmp = a[0]; a[0] = b[0]; b[0] = tmp; \
tmp = a[1]; a[1] = b[1]; b[1] = tmp; \
tmp = a[2]; a[2] = b[2]; b[2] = tmp; \
} while(0)
/* ---------------------------------------------------------------------- */
void ComputeOrientOrderAtom::select3(int k, int n, double *arr, int *iarr, double **arr3)
{
int i,ir,j,l,mid,ia,itmp;
double a,tmp,a3[3];
arr--;
iarr--;
arr3--;
l = 1;
ir = n;
for (;;) {
if (ir <= l+1) {
if (ir == l+1 && arr[ir] < arr[l]) {
SWAP(arr[l],arr[ir]);
ISWAP(iarr[l],iarr[ir]);
SWAP3(arr3[l],arr3[ir]);
}
return;
} else {
mid=(l+ir) >> 1;
SWAP(arr[mid],arr[l+1]);
ISWAP(iarr[mid],iarr[l+1]);
SWAP3(arr3[mid],arr3[l+1]);
if (arr[l] > arr[ir]) {
SWAP(arr[l],arr[ir]);
ISWAP(iarr[l],iarr[ir]);
SWAP3(arr3[l],arr3[ir]);
}
if (arr[l+1] > arr[ir]) {
SWAP(arr[l+1],arr[ir]);
ISWAP(iarr[l+1],iarr[ir]);
SWAP3(arr3[l+1],arr3[ir]);
}
if (arr[l] > arr[l+1]) {
SWAP(arr[l],arr[l+1]);
ISWAP(iarr[l],iarr[l+1]);
SWAP3(arr3[l],arr3[l+1]);
}
i = l+1;
j = ir;
a = arr[l+1];
ia = iarr[l+1];
a3[0] = arr3[l+1][0];
a3[1] = arr3[l+1][1];
a3[2] = arr3[l+1][2];
for (;;) {
do i++; while (arr[i] < a);
do j--; while (arr[j] > a);
if (j < i) break;
SWAP(arr[i],arr[j]);
ISWAP(iarr[i],iarr[j]);
SWAP3(arr3[i],arr3[j]);
}
arr[l+1] = arr[j];
arr[j] = a;
iarr[l+1] = iarr[j];
iarr[j] = ia;
arr3[l+1][0] = arr3[j][0];
arr3[l+1][1] = arr3[j][1];
arr3[l+1][2] = arr3[j][2];
arr3[j][0] = a3[0];
arr3[j][1] = a3[1];
arr3[j][2] = a3[2];
if (j >= k) ir = j-1;
if (j <= k) l = i;
}
}
}
/* ----------------------------------------------------------------------
calculate the bond orientational order parameters
------------------------------------------------------------------------- */
void ComputeOrientOrderAtom::calc_boop(double **rlist,
int numNeighbors, double qn[],
int qlist[], int nqlist) {
// DONE: Need to make this memory allocation dynamic,
// DONE: Add error handling
// DONE: Add to memory usage
// DONE: Add options for different nnn and qlist
for (int iw = 0; iw < nqlist; iw++) {
int n = qlist[iw];
qn[iw] = 0.0;
for(int m = 0; m < 2*n+1; m++) {
qnm_r[iw][m] = 0.0;
qnm_i[iw][m] = 0.0;
}
}
for(int ineigh = 0; ineigh < numNeighbors; ineigh++) {
const double * const r = rlist[ineigh];
double rmag = dist(r);
if(rmag <= MY_EPSILON) {
return;
}
double costheta = r[2] / rmag;
double expphi_r = r[0];
double expphi_i = r[1];
double rxymag = sqrt(expphi_r*expphi_r+expphi_i*expphi_i);
if(rxymag <= MY_EPSILON) {
expphi_r = 1.0;
expphi_i = 0.0;
} else {
double rxymaginv = 1.0/rxymag;
expphi_r *= rxymaginv;
expphi_i *= rxymaginv;
}
for (int iw = 0; iw < nqlist; iw++) {
int n = qlist[iw];
qnm_r[iw][n] += polar_prefactor(n, 0, costheta);
double expphim_r = expphi_r;
double expphim_i = expphi_i;
for(int m = 1; m <= +n; m++) {
double prefactor = polar_prefactor(n, m, costheta);
double c_r = prefactor * expphim_r;
double c_i = prefactor * expphim_i;
qnm_r[iw][m+n] += c_r;
qnm_i[iw][m+n] += c_i;
if(m & 1) {
qnm_r[iw][-m+n] -= c_r;
qnm_i[iw][-m+n] += c_i;
} else {
qnm_r[iw][-m+n] += c_r;
qnm_i[iw][-m+n] -= c_i;
}
double tmp_r = expphim_r*expphi_r - expphim_i*expphi_i;
double tmp_i = expphim_r*expphi_i + expphim_i*expphi_r;
expphim_r = tmp_r;
expphim_i = tmp_i;
}
}
}
double fac = sqrt(MY_4PI) / numNeighbors;
for (int iw = 0; iw < nqlist; iw++) {
int n = qlist[iw];
double qm_sum = 0.0;
for(int m = 0; m < 2*n+1; m++) {
qm_sum += qnm_r[iw][m]*qnm_r[iw][m] + qnm_i[iw][m]*qnm_i[iw][m];
}
qn[iw] = fac * sqrt(qm_sum / (2*n+1));
}
}
/* ----------------------------------------------------------------------
calculate scalar distance
------------------------------------------------------------------------- */
double ComputeOrientOrderAtom::dist(const double r[]) {
return sqrt(r[0]*r[0] + r[1]*r[1] + r[2]*r[2]);
}
/* ----------------------------------------------------------------------
polar prefactor for spherical harmonic Y_l^m, where
Y_l^m (theta, phi) = prefactor(l, m, cos(theta)) * exp(i*m*phi)
------------------------------------------------------------------------- */
double ComputeOrientOrderAtom::
polar_prefactor(int l, int m, double costheta) {
const int mabs = abs(m);
double prefactor = 1.0;
for (int i=l-mabs+1; i < l+mabs+1; ++i)
prefactor *= static_cast<double>(i);
prefactor = sqrt(static_cast<double>(2*l+1)/(MY_4PI*prefactor))
* associated_legendre(l,mabs,costheta);
if ((m < 0) && (m % 2)) prefactor = -prefactor;
return prefactor;
}
/* ----------------------------------------------------------------------
associated legendre polynomial
------------------------------------------------------------------------- */
double ComputeOrientOrderAtom::
associated_legendre(int l, int m, double x) {
if (l < m) return 0.0;
double p(1.0), pm1(0.0), pm2(0.0);
if (m != 0) {
const double sqx = sqrt(1.0-x*x);
for (int i=1; i < m+1; ++i)
p *= static_cast<double>(2*i-1) * sqx;
}
for (int i=m+1; i < l+1; ++i) {
pm2 = pm1;
pm1 = p;
p = (static_cast<double>(2*i-1)*x*pm1
- static_cast<double>(i+m-1)*pm2) / static_cast<double>(i-m);
}
return p;
}

View File

@ -0,0 +1,80 @@
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#ifdef COMPUTE_CLASS
ComputeStyle(orientorder/atom,ComputeOrientOrderAtom)
#else
#ifndef LMP_COMPUTE_ORIENTORDER_ATOM_H
#define LMP_COMPUTE_ORIENTORDER_ATOM_H
#include "compute.h"
namespace LAMMPS_NS {
class ComputeOrientOrderAtom : public Compute {
public:
ComputeOrientOrderAtom(class LAMMPS *, int, char **);
~ComputeOrientOrderAtom();
void init();
void init_list(int, class NeighList *);
void compute_peratom();
double memory_usage();
private:
int nmax,maxneigh,ncol,nnn;
class NeighList *list;
double *distsq;
int *nearest;
double **rlist;
int *qlist;
int nqlist;
int qmax;
double **qnarray;
double cutsq;
double **qnm_r;
double **qnm_i;
void select3(int, int, double *, int *, double **);
void calc_boop(double **rlist, int numNeighbors,
double qn[], int nlist[], int nnlist);
double dist(const double r[]);
double polar_prefactor(int, int, double);
double associated_legendre(int, int, double);
};
}
#endif
#endif
/* ERROR/WARNING messages:
E: Illegal ... command
Self-explanatory. Check the input script syntax and compare to the
documentation for the command. You can use -echo screen as a
command-line option when running LAMMPS to see the offending line.
E: Compute orientorder/atom requires a pair style be defined
Self-explantory.
W: More than one compute orientorder/atom
It is not efficient to use compute orientorder/atom more than once.
*/