forked from lijiext/lammps
415 lines
19 KiB
HTML
415 lines
19 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>pair_style hbond/dreiding/lj command — 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 & scalability</a></li>
|
|
<li class="toctree-l1"><a class="reference internal" href="Section_tools.html">9. Additional tools</a></li>
|
|
<li class="toctree-l1"><a class="reference internal" href="Section_modify.html">10. Modifying & extending LAMMPS</a></li>
|
|
<li class="toctree-l1"><a class="reference internal" href="Section_python.html">11. Python interface to LAMMPS</a></li>
|
|
<li class="toctree-l1"><a class="reference internal" href="Section_errors.html">12. Errors</a></li>
|
|
<li class="toctree-l1"><a class="reference internal" href="Section_history.html">13. Future and history</a></li>
|
|
</ul>
|
|
|
|
|
|
|
|
</div>
|
|
|
|
</nav>
|
|
|
|
<section data-toggle="wy-nav-shift" class="wy-nav-content-wrap">
|
|
|
|
|
|
<nav class="wy-nav-top" role="navigation" aria-label="top navigation">
|
|
<i data-toggle="wy-nav-top" class="fa fa-bars"></i>
|
|
<a href="Manual.html">LAMMPS</a>
|
|
</nav>
|
|
|
|
|
|
|
|
<div class="wy-nav-content">
|
|
<div class="rst-content">
|
|
<div role="navigation" aria-label="breadcrumbs navigation">
|
|
<ul class="wy-breadcrumbs">
|
|
<li><a href="Manual.html">Docs</a> »</li>
|
|
|
|
<li>pair_style hbond/dreiding/lj 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="pair-style-hbond-dreiding-lj-command">
|
|
<span id="index-0"></span><h1>pair_style hbond/dreiding/lj command<a class="headerlink" href="#pair-style-hbond-dreiding-lj-command" title="Permalink to this headline">¶</a></h1>
|
|
</div>
|
|
<div class="section" id="pair-style-hbond-dreiding-lj-omp-command">
|
|
<h1>pair_style hbond/dreiding/lj/omp command<a class="headerlink" href="#pair-style-hbond-dreiding-lj-omp-command" title="Permalink to this headline">¶</a></h1>
|
|
</div>
|
|
<div class="section" id="pair-style-hbond-dreiding-morse-command">
|
|
<h1>pair_style hbond/dreiding/morse command<a class="headerlink" href="#pair-style-hbond-dreiding-morse-command" title="Permalink to this headline">¶</a></h1>
|
|
</div>
|
|
<div class="section" id="pair-style-hbond-dreiding-morse-omp-command">
|
|
<h1>pair_style hbond/dreiding/morse/omp command<a class="headerlink" href="#pair-style-hbond-dreiding-morse-omp-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>pair_style style N inner_distance_cutoff outer_distance_cutoff angle_cutof
|
|
</pre></div>
|
|
</div>
|
|
<ul class="simple">
|
|
<li>style = <em>hbond/dreiding/lj</em> or <em>hbond/dreiding/morse</em></li>
|
|
<li>n = cosine angle periodicity</li>
|
|
<li>inner_distance_cutoff = global inner cutoff for Donor-Acceptor interactions (distance units)</li>
|
|
<li>outer_distance_cutoff = global cutoff for Donor-Acceptor interactions (distance units)</li>
|
|
<li>angle_cutoff = global angle cutoff for Acceptor-Hydrogen-Donor</li>
|
|
<li>interactions (degrees)</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>pair_style hybrid/overlay lj/cut 10.0 hbond/dreiding/lj 4 9.0 11.0 90
|
|
pair_coeff 1 2 hbond/dreiding/lj 3 i 9.5 2.75 4 9.0 11.0 90.0
|
|
</pre></div>
|
|
</div>
|
|
<div class="highlight-python"><div class="highlight"><pre>pair_style hybrid/overlay lj/cut 10.0 hbond/dreiding/morse 2 9.0 11.0 90
|
|
pair_coeff 1 2 hbond/dreiding/morse 3 i 3.88 1.7241379 2.9 2 9 11 90
|
|
</pre></div>
|
|
</div>
|
|
</div>
|
|
<div class="section" id="description">
|
|
<h2>Description<a class="headerlink" href="#description" title="Permalink to this headline">¶</a></h2>
|
|
<p>The <em>hbond/dreiding</em> styles compute the Acceptor-Hydrogen-Donor (AHD)
|
|
3-body hydrogen bond interaction for the
|
|
<a class="reference internal" href="Section_howto.html#howto-4"><span>DREIDING</span></a> force field, given by:</p>
|
|
<img alt="_images/pair_hbond_dreiding.jpg" class="align-center" src="_images/pair_hbond_dreiding.jpg" />
|
|
<p>where Rin is the inner spline distance cutoff, Rout is the outer
|
|
distance cutoff, theta_c is the angle cutoff, and n is the cosine
|
|
periodicity.</p>
|
|
<p>Here, <em>r</em> is the radial distance between the donor (D) and acceptor
|
|
(A) atoms and <em>theta</em> is the bond angle between the acceptor, the
|
|
hydrogen (H) and the donor atoms:</p>
|
|
<img alt="_images/dreiding_hbond.jpg" class="align-center" src="_images/dreiding_hbond.jpg" />
|
|
<p>These 3-body interactions can be defined for pairs of acceptor and
|
|
donor atoms, based on atom types. For each donor/acceptor atom pair,
|
|
the 3rd atom in the interaction is a hydrogen permanently bonded to
|
|
the donor atom, e.g. in a bond list read in from a data file via the
|
|
<a class="reference internal" href="read_data.html"><em>read_data</em></a> command. The atom types of possible
|
|
hydrogen atoms for each donor/acceptor type pair are specified by the
|
|
<a class="reference internal" href="pair_coeff.html"><em>pair_coeff</em></a> command (see below).</p>
|
|
<p>Style <em>hbond/dreiding/lj</em> is the original DREIDING potential of
|
|
<a class="reference internal" href="special_bonds.html#mayo"><span>(Mayo)</span></a>. It uses a LJ 12/10 functional for the Donor-Acceptor
|
|
interactions. To match the results in the original paper, use n = 4.</p>
|
|
<p>Style <em>hbond/dreiding/morse</em> is an improved version using a Morse
|
|
potential for the Donor-Acceptor interactions. <a class="reference internal" href="#liu"><span>(Liu)</span></a> showed
|
|
that the Morse form gives improved results for Dendrimer simulations,
|
|
when n = 2.</p>
|
|
<p>See this <a class="reference internal" href="Section_howto.html#howto-4"><span>howto section</span></a> of the manual for
|
|
more information on the DREIDING forcefield.</p>
|
|
<div class="admonition warning">
|
|
<p class="first admonition-title">Warning</p>
|
|
<p class="last">Because the Dreiding hydrogen bond potential is only
|
|
one portion of an overall force field which typically includes other
|
|
pairwise interactions, it is common to use it as a sub-style in a
|
|
<a class="reference internal" href="pair_hybrid.html"><em>pair_style hybrid/overlay</em></a> command, where another
|
|
pair style provides the repulsive core interaction between pairs of
|
|
atoms, e.g. a 1/r^12 Lennard-Jones repulsion.</p>
|
|
</div>
|
|
<div class="admonition warning">
|
|
<p class="first admonition-title">Warning</p>
|
|
<p class="last">When using the hbond/dreiding pair styles with
|
|
<a class="reference internal" href="pair_hybrid.html"><em>pair_style hybrid/overlay</em></a>, you should explicitly
|
|
define pair interactions between the donor atom and acceptor atoms,
|
|
(as well as between these atoms and ALL other atoms in your system).
|
|
Whenever <a class="reference internal" href="pair_hybrid.html"><em>pair_style hybrid/overlay</em></a> is used,
|
|
ordinary mixing rules are not applied to atoms like the donor and
|
|
acceptor atoms because they are typically referenced in multiple pair
|
|
styles. Neglecting to do this can cause difficult-to-detect physics
|
|
problems.</p>
|
|
</div>
|
|
<div class="admonition warning">
|
|
<p class="first admonition-title">Warning</p>
|
|
<p class="last">In the original Dreiding force field paper 1-4
|
|
non-bonded interactions ARE allowed. If this is desired for your
|
|
model, use the special_bonds command (e.g. “special_bonds lj 0.0 0.0
|
|
1.0”) to turn these interactions on.</p>
|
|
</div>
|
|
<hr class="docutils" />
|
|
<p>The following coefficients must be defined for pairs of eligible
|
|
donor/acceptor types via the <a class="reference internal" href="pair_coeff.html"><em>pair_coeff</em></a> command as
|
|
in the examples above.</p>
|
|
<div class="admonition warning">
|
|
<p class="first admonition-title">Warning</p>
|
|
<p class="last">Unlike other pair styles and their associated
|
|
<a class="reference internal" href="pair_coeff.html"><em>pair_coeff</em></a> commands, you do not need to specify
|
|
pair_coeff settings for all possible I,J type pairs. Only I,J type
|
|
pairs for atoms which act as joint donors/acceptors need to be
|
|
specified; all other type pairs are assumed to be inactive.</p>
|
|
</div>
|
|
<div class="admonition warning">
|
|
<p class="first admonition-title">Warning</p>
|
|
<p class="last">A <a class="reference internal" href="pair_coeff.html"><em>pair_coeff</em></a> command can be
|
|
speficied multiple times for the same donor/acceptor type pair. This
|
|
enables multiple hydrogen types to be assigned to the same
|
|
donor/acceptor type pair. For other pair_styles, if the pair_coeff
|
|
command is re-used for the same I.J type pair, the settings for that
|
|
type pair are overwritten. For the hydrogen bond potentials this is
|
|
not the case; the settings are cummulative. This means the only way
|
|
to turn off a previous setting, is to re-use the pair_style command
|
|
and start over.</p>
|
|
</div>
|
|
<p>For the <em>hbond/dreiding/lj</em> style the list of coefficients is as
|
|
follows:</p>
|
|
<ul class="simple">
|
|
<li>K = hydrogen atom type = 1 to Ntypes</li>
|
|
<li>donor flag = <em>i</em> or <em>j</em></li>
|
|
<li>epsilon (energy units)</li>
|
|
<li>sigma (distance units)</li>
|
|
<li>n = exponent in formula above</li>
|
|
<li>distance cutoff Rin (distance units)</li>
|
|
<li>distance cutoff Rout (distance units)</li>
|
|
<li>angle cutoff (degrees)</li>
|
|
</ul>
|
|
<p>For the <em>hbond/dreiding/morse</em> style the list of coefficients is as
|
|
follows:</p>
|
|
<ul class="simple">
|
|
<li>K = hydrogen atom type = 1 to Ntypes</li>
|
|
<li>donor flag = <em>i</em> or <em>j</em></li>
|
|
<li>D0 (energy units)</li>
|
|
<li>alpha (1/distance units)</li>
|
|
<li>r0 (distance units)</li>
|
|
<li>n = exponent in formula above</li>
|
|
<li>distance cutoff Rin (distance units)</li>
|
|
<li>distance cutoff Rout (distance units)</li>
|
|
<li>angle cutoff (degrees)</li>
|
|
</ul>
|
|
<p>A single hydrogen atom type K can be specified, or a wild-card
|
|
asterisk can be used in place of or in conjunction with the K
|
|
arguments to select multiple types as hydrogens. This takes the form
|
|
“*” or “<em>n” or “n</em>” or “m*n”. See the <a class="reference external" href="pair_coeff">pair_coeff</a> command
|
|
doc page for details.</p>
|
|
<p>If the donor flag is <em>i</em>, then the atom of type I in the pair_coeff
|
|
command is treated as the donor, and J is the acceptor. If the donor
|
|
flag is <em>j</em>, then the atom of type J in the pair_coeff command is
|
|
treated as the donor and I is the donor. This option is required
|
|
because the <a class="reference internal" href="pair_coeff.html"><em>pair_coeff</em></a> command requires that I <= J.</p>
|
|
<p>Epsilon and sigma are settings for the hydrogen bond potential based
|
|
on a Lennard-Jones functional form. Note that sigma is defined as the
|
|
zero-crossing distance for the potential, not as the energy minimum at
|
|
2^(1/6) sigma.</p>
|
|
<p>D0 and alpha and r0 are settings for the hydrogen bond potential based
|
|
on a Morse functional form.</p>
|
|
<p>The last 3 coefficients for both styles are optional. If not
|
|
specified, the global n, distance cutoff, and angle cutoff specified
|
|
in the pair_style command are used. If you wish to only override the
|
|
2nd or 3rd optional parameter, you must also specify the preceding
|
|
optional parameters.</p>
|
|
<hr class="docutils" />
|
|
<p>Styles with a <em>cuda</em>, <em>gpu</em>, <em>intel</em>, <em>kk</em>, <em>omp</em>, or <em>opt</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, GPU, USER-INTEL,
|
|
KOKKOS, USER-OMP and OPT packages, respectively. They are only
|
|
enabled if LAMMPS was built with those packages. 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>
|
|
<hr class="docutils" />
|
|
<p><strong>Mixing, shift, table, tail correction, restart, rRESPA info</strong>:</p>
|
|
<p>These pair styles do not support mixing. You must explicitly identify
|
|
each donor/acceptor type pair.</p>
|
|
<p>These styles do not support the <a class="reference internal" href="pair_modify.html"><em>pair_modify</em></a> shift
|
|
option for the energy of the interactions.</p>
|
|
<p>The <a class="reference internal" href="pair_modify.html"><em>pair_modify</em></a> table option is not relevant for
|
|
these pair styles.</p>
|
|
<p>These pair styles do not support the <a class="reference internal" href="pair_modify.html"><em>pair_modify</em></a>
|
|
tail option for adding long-range tail corrections to energy and
|
|
pressure.</p>
|
|
<p>These pair styles do not write their information to <a class="reference internal" href="restart.html"><em>binary restart files</em></a>, so pair_style and pair_coeff commands need to be
|
|
re-specified in an input script that reads a restart file.</p>
|
|
<p>These pair styles can only be used via the <em>pair</em> keyword of the
|
|
<a class="reference internal" href="run_style.html"><em>run_style respa</em></a> command. They do not support the
|
|
<em>inner</em>, <em>middle</em>, <em>outer</em> keywords.</p>
|
|
<p>These pair styles tally a count of how many hydrogen bonding
|
|
interactions they calculate each timestep and the hbond energy. These
|
|
quantities can be accessed via the <a class="reference internal" href="compute_pair.html"><em>compute pair</em></a>
|
|
command as a vector of values of length 2.</p>
|
|
<p>To print these quantities to the log file (with a descriptive column
|
|
heading) the following commands could be included in an input script:</p>
|
|
<div class="highlight-python"><div class="highlight"><pre>compute hb all pair hbond/dreiding/lj
|
|
variable n_hbond equal c_hb[1] #number hbonds
|
|
variable E_hbond equal c_hb[2] #hbond energy
|
|
thermo_style custom step temp epair v_E_hbond
|
|
</pre></div>
|
|
</div>
|
|
</div>
|
|
<hr class="docutils" />
|
|
<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="pair_coeff.html"><em>pair_coeff</em></a></p>
|
|
<p><strong>Default:</strong> none</p>
|
|
<hr class="docutils" />
|
|
<p id="mayo"><strong>(Mayo)</strong> Mayo, Olfason, Goddard III, J Phys Chem, 94, 8897-8909
|
|
(1990).</p>
|
|
<p id="liu"><strong>(Liu)</strong> Liu, Bryantsev, Diallo, Goddard III, J. Am. Chem. Soc 131 (8)
|
|
2798 (2009)</p>
|
|
</div>
|
|
</div>
|
|
|
|
|
|
</div>
|
|
</div>
|
|
<footer>
|
|
|
|
|
|
<hr/>
|
|
|
|
<div role="contentinfo">
|
|
<p>
|
|
© 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> |