2009-11-07 05:01:28 +08:00
2015-07-30 22:53:28 +08:00
<!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 yukawa/colloid command — LAMMPS 15 May 2015 version documentation< / title >
2009-11-07 05:01:28 +08:00
2015-07-30 22:53:28 +08:00
2009-11-07 05:01:28 +08:00
2015-07-30 22:53:28 +08:00
2009-11-07 05:01:28 +08:00
2015-07-30 22:53:28 +08:00
2009-11-07 05:01:28 +08:00
2015-07-30 22:53:28 +08:00
2009-11-07 05:01:28 +08:00
2015-07-30 22:53:28 +08:00
< 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 yukawa/colloid 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-yukawa-colloid-command" >
< span id = "index-0" > < / span > < h1 > pair_style yukawa/colloid command< a class = "headerlink" href = "#pair-style-yukawa-colloid-command" title = "Permalink to this headline" > ¶< / a > < / h1 >
< / div >
< div class = "section" id = "pair-style-yukawa-colloid-gpu-command" >
< h1 > pair_style yukawa/colloid/gpu command< a class = "headerlink" href = "#pair-style-yukawa-colloid-gpu-command" title = "Permalink to this headline" > ¶< / a > < / h1 >
< / div >
< div class = "section" id = "pair-style-yukawa-colloid-omp-command" >
< h1 > pair_style yukawa/colloid/omp command< a class = "headerlink" href = "#pair-style-yukawa-colloid-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 yukawa/colloid kappa cutoff
< / pre > < / div >
< / div >
< ul class = "simple" >
< li > kappa = screening length (inverse distance units)< / li >
< li > cutoff = global cutoff for colloidal Yukawa interactions (distance units)< / li >
< / ul >
< / div >
< div class = "section" id = "examples" >
< h2 > Examples< a class = "headerlink" href = "#examples" title = "Permalink to this headline" > ¶< / a > < / h2 >
< div class = "highlight-python" > < div class = "highlight" > < pre > pair_style yukawa/colloid 2.0 2.5
2009-11-07 05:01:28 +08:00
pair_coeff 1 1 100.0 2.3
2015-07-30 22:53:28 +08:00
pair_coeff * * 100.0
< / pre > < / div >
< / div >
< / div >
< div class = "section" id = "description" >
< h2 > Description< a class = "headerlink" href = "#description" title = "Permalink to this headline" > ¶< / a > < / h2 >
< p > Style < em > yukawa/colloid< / em > computes pairwise interactions with the formula< / p >
< img alt = "_images/pair_yukawa_colloid.jpg" class = "align-center" src = "_images/pair_yukawa_colloid.jpg" / >
< p > where Ri and Rj are the radii of the two particles and Rc is the
cutoff.< / p >
< p > In contrast to < a class = "reference internal" href = "pair_yukawa.html" > < em > pair_style yukawa< / em > < / a > , this functional
2009-11-07 05:01:28 +08:00
form arises from the Coulombic interaction between two colloid
2013-07-02 22:18:40 +08:00
particles, screened due to the presence of an electrolyte, see the
2015-07-30 22:53:28 +08:00
book by < a class = "reference internal" href = "#safran" > < span > Safran< / span > < / a > for a derivation in the context of DVLO
theory. < a class = "reference internal" href = "pair_yukawa.html" > < em > Pair_style yukawa< / em > < / a > is a screened Coulombic
potential between two point-charges and uses no such approximation.< / p >
< p > This potential applies to nearby particle pairs for which the Derjagin
approximation holds, meaning h < < Ri + Rj, where h is the
surface-to-surface separation of the two particles.< / p >
< p > When used in combination with < a class = "reference internal" href = "pair_colloid.html" > < em > pair_style colloid< / em > < / a > ,
2009-11-07 05:01:28 +08:00
the two terms become the so-called DLVO potential, which combines
2015-07-30 22:53:28 +08:00
electrostatic repulsion and van der Waals attraction.< / p >
< p > The following coefficients must be defined for each pair of atoms
types via the < a class = "reference internal" href = "pair_coeff.html" > < em > pair_coeff< / em > < / a > command as in the examples
2009-11-07 05:01:28 +08:00
above, or in the data file or restart files read by the
2015-07-30 22:53:28 +08:00
< a class = "reference internal" href = "read_data.html" > < em > read_data< / em > < / a > or < a class = "reference internal" href = "read_restart.html" > < em > read_restart< / em > < / a >
commands, or by mixing as described below:< / p >
< ul class = "simple" >
< li > A (energy/distance units)< / li >
< li > cutoff (distance units)< / li >
< / ul >
< p > The prefactor A is determined from the relationship between surface
2009-11-07 05:01:28 +08:00
charge and surface potential due to the presence of electrolyte. Note
that the A for this potential style has different units than the A
2015-07-30 22:53:28 +08:00
used in < a class = "reference internal" href = "pair_yukawa.html" > < em > pair_style yukawa< / em > < / a > . For low surface
potentials, i.e. less than about 25 mV, A can be written as:< / p >
< div class = "highlight-python" > < div class = "highlight" > < pre > < span class = "n" > A< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "o" > *< / span > < span class = "n" > PI< / span > < span class = "o" > *< / span > < span class = "n" > R< / span > < span class = "o" > *< / span > < span class = "n" > eps< / span > < span class = "o" > *< / span > < span class = "n" > eps0< / span > < span class = "o" > *< / span > < span class = "n" > kappa< / span > < span class = "o" > *< / span > < span class = "n" > psi< / span > < span class = "o" > ^< / span > < span class = "mi" > 2< / span >
< / pre > < / div >
< / div >
< p > where< / p >
< ul class = "simple" >
< li > R = colloid radius (distance units)< / li >
< li > eps0 = permittivity of free space (charge^2/energy/distance units)< / li >
< li > eps = relative permittivity of fluid medium (dimensionless)< / li >
< li > kappa = inverse screening length (1/distance units)< / li >
< li > psi = surface potential (energy/charge units)< / li >
< / ul >
< p > The last coefficient is optional. If not specified, the global
yukawa/colloid cutoff is used.< / 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
2014-08-15 00:30:25 +08:00
functionally the same as the corresponding style without the suffix.
They have been optimized to run faster, depending on your available
2015-07-30 22:53:28 +08:00
hardware, as discussed in < a class = "reference internal" href = "Section_accelerate.html" > < em > Section_accelerate< / em > < / a >
2014-08-15 00:30:25 +08:00
of the manual. The accelerated styles take the same arguments and
should produce the same results, except for round-off and precision
2015-07-30 22:53:28 +08:00
issues.< / p >
< p > These accelerated styles are part of the USER-CUDA, GPU, USER-INTEL,
2014-08-15 00:30:25 +08:00
KOKKOS, USER-OMP and OPT packages, respectively. They are only
2015-07-30 22:53:28 +08:00
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 > For atom type pairs I,J and I != J, the A coefficient and cutoff
2009-11-07 05:01:28 +08:00
distance for this pair style can be mixed. A is an energy value mixed
2015-07-30 22:53:28 +08:00
like a LJ epsilon. The default mix value is < em > geometric< / em > . See the
“ pair_modify” command for details.< / p >
< p > This pair style supports the < a class = "reference internal" href = "pair_modify.html" > < em > pair_modify< / em > < / a > shift
option for the energy of the pair interaction.< / p >
< p > The < a class = "reference internal" href = "pair_modify.html" > < em > pair_modify< / em > < / a > table option is not relevant
for this pair style.< / p >
< p > This pair style does not support the < a class = "reference internal" href = "pair_modify.html" > < em > pair_modify< / em > < / a >
2009-11-07 05:01:28 +08:00
tail option for adding long-range tail corrections to energy and
2015-07-30 22:53:28 +08:00
pressure.< / p >
< p > This pair style writes its information to < a class = "reference internal" href = "restart.html" > < em > binary restart files< / em > < / a > , so pair_style and pair_coeff commands do not need
to be specified in an input script that reads a restart file.< / p >
< p > This pair style 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. It does not support the
< em > inner< / em > , < em > middle< / em > , < em > outer< / em > keywords.< / p >
< / div >
< hr class = "docutils" / >
< div class = "section" id = "restrictions" >
< h2 > Restrictions< a class = "headerlink" href = "#restrictions" title = "Permalink to this headline" > ¶< / a > < / h2 >
< p > This style is part of the COLLOID 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 > Making LAMMPS< / span > < / a > section for more info.< / p >
< p > This pair style requires that atoms be finite-size spheres with a
diameter, as defined by the < a class = "reference internal" href = "atom_style.html" > < em > atom_style sphere< / em > < / a >
command.< / p >
< p > Per-particle polydispersity is not yet supported by this pair style;
2011-04-14 05:39:34 +08:00
per-type polydispersity is allowed. This means all particles of the
same type must have the same diameter. Each type can have a different
2015-07-30 22:53:28 +08:00
diameter.< / p >
< / div >
< div class = "section" id = "related-commands" >
< h2 > Related commands< a class = "headerlink" href = "#related-commands" title = "Permalink to this headline" > ¶< / a > < / h2 >
< p > < a class = "reference internal" href = "pair_coeff.html" > < em > pair_coeff< / em > < / a > < / p >
< p > < strong > Default:< / strong > none< / p >
< hr class = "docutils" / >
< p id = "safran" > < strong > (Safran)< / strong > Safran, Statistical Thermodynamics of Surfaces, Interfaces,
And Membranes, Westview Press, ISBN: 978-0813340791 (2003).< / 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 >