Files
lammps/doc/html/Howto_coreshell.html
2025-01-13 14:55:48 +00:00

393 lines
34 KiB
HTML

<!DOCTYPE html>
<html class="writer-html5" lang="en" >
<head>
<meta charset="utf-8" /><meta name="viewport" content="width=device-width, initial-scale=1" />
<meta name="viewport" content="width=device-width, initial-scale=1.0" />
<title>8.5.6. Adiabatic core/shell model &mdash; LAMMPS documentation</title>
<link rel="stylesheet" href="_static/pygments.css" type="text/css" />
<link rel="stylesheet" href="_static/css/theme.css" type="text/css" />
<link rel="stylesheet" href="_static/sphinx-design.min.css" type="text/css" />
<link rel="stylesheet" href="_static/css/lammps.css" type="text/css" />
<link rel="shortcut icon" href="_static/lammps.ico"/>
<link rel="canonical" href="https://docs.lammps.org/Howto_coreshell.html" />
<!--[if lt IE 9]>
<script src="_static/js/html5shiv.min.js"></script>
<![endif]-->
<script src="_static/jquery.js?v=5d32c60e"></script>
<script src="_static/_sphinx_javascript_frameworks_compat.js?v=2cd50e6c"></script>
<script src="_static/documentation_options.js?v=5929fcd5"></script>
<script src="_static/doctools.js?v=9bcbadda"></script>
<script src="_static/sphinx_highlight.js?v=dc90522c"></script>
<script src="_static/design-tabs.js?v=f930bc37"></script>
<script async="async" src="_static/mathjax/es5/tex-mml-chtml.js?v=cadf963e"></script>
<script src="_static/js/theme.js"></script>
<link rel="index" title="Index" href="genindex.html" />
<link rel="search" title="Search" href="search.html" />
<link rel="next" title="8.5.7. Drude induced dipoles" href="Howto_drude.html" />
<link rel="prev" title="8.5.5. Polarizable models" href="Howto_polarizable.html" />
</head>
<body class="wy-body-for-nav">
<div class="wy-grid-for-nav">
<nav data-toggle="wy-nav-shift" class="wy-nav-side">
<div class="wy-side-scroll">
<div class="wy-side-nav-search" >
<a href="Manual.html">
<img src="_static/lammps-logo.png" class="logo" alt="Logo"/>
</a>
<div class="lammps_version">Version: <b>19 Nov 2024</b></div>
<div class="lammps_release">git info: </div>
<div role="search">
<form id="rtd-search-form" class="wy-form" action="search.html" method="get">
<input type="text" name="q" placeholder="Search docs" aria-label="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="Navigation menu">
<p class="caption" role="heading"><span class="caption-text">User Guide</span></p>
<ul class="current">
<li class="toctree-l1"><a class="reference internal" href="Intro.html">1. Introduction</a></li>
<li class="toctree-l1"><a class="reference internal" href="Install.html">2. Install LAMMPS</a></li>
<li class="toctree-l1"><a class="reference internal" href="Build.html">3. Build LAMMPS</a></li>
<li class="toctree-l1"><a class="reference internal" href="Run_head.html">4. Run LAMMPS</a></li>
<li class="toctree-l1"><a class="reference internal" href="Commands.html">5. Commands</a></li>
<li class="toctree-l1"><a class="reference internal" href="Packages.html">6. Optional packages</a></li>
<li class="toctree-l1"><a class="reference internal" href="Speed.html">7. Accelerate performance</a></li>
<li class="toctree-l1 current"><a class="reference internal" href="Howto.html">8. Howto discussions</a><ul class="current">
<li class="toctree-l2"><a class="reference internal" href="Howto.html#general-howto">8.1. General howto</a></li>
<li class="toctree-l2"><a class="reference internal" href="Howto.html#settings-howto">8.2. Settings howto</a></li>
<li class="toctree-l2"><a class="reference internal" href="Howto.html#analysis-howto">8.3. Analysis howto</a></li>
<li class="toctree-l2"><a class="reference internal" href="Howto.html#force-fields-howto">8.4. Force fields howto</a></li>
<li class="toctree-l2 current"><a class="reference internal" href="Howto.html#packages-howto">8.5. Packages howto</a><ul class="current">
<li class="toctree-l3"><a class="reference internal" href="Howto_spherical.html">8.5.1. Finite-size spherical and aspherical particles</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_granular.html">8.5.2. Granular models</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_body.html">8.5.3. Body particles</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_bpm.html">8.5.4. Bonded particle models</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_polarizable.html">8.5.5. Polarizable models</a></li>
<li class="toctree-l3 current"><a class="current reference internal" href="#">8.5.6. Adiabatic core/shell model</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_drude.html">8.5.7. Drude induced dipoles</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_drude2.html">8.5.8. Tutorial for Thermalized Drude oscillators in LAMMPS</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_peri.html">8.5.9. Peridynamics with LAMMPS</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_manifold.html">8.5.10. Manifolds (surfaces)</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_rheo.html">8.5.11. Reproducing hydrodynamics and elastic objects (RHEO)</a></li>
<li class="toctree-l3"><a class="reference internal" href="Howto_spins.html">8.5.12. Magnetic spins</a></li>
</ul>
</li>
<li class="toctree-l2"><a class="reference internal" href="Howto.html#tutorials-howto">8.6. Tutorials howto</a></li>
</ul>
</li>
<li class="toctree-l1"><a class="reference internal" href="Examples.html">9. Example scripts</a></li>
<li class="toctree-l1"><a class="reference internal" href="Tools.html">10. Auxiliary tools</a></li>
<li class="toctree-l1"><a class="reference internal" href="Errors.html">11. Errors</a></li>
</ul>
<p class="caption" role="heading"><span class="caption-text">Programmer Guide</span></p>
<ul>
<li class="toctree-l1"><a class="reference internal" href="Library.html">1. LAMMPS Library Interfaces</a></li>
<li class="toctree-l1"><a class="reference internal" href="Python_head.html">2. Use Python with LAMMPS</a></li>
<li class="toctree-l1"><a class="reference internal" href="Modify.html">3. Modifying &amp; extending LAMMPS</a></li>
<li class="toctree-l1"><a class="reference internal" href="Developer.html">4. Information for Developers</a></li>
</ul>
<p class="caption" role="heading"><span class="caption-text">Command Reference</span></p>
<ul>
<li class="toctree-l1"><a class="reference internal" href="commands_list.html">Commands</a></li>
<li class="toctree-l1"><a class="reference internal" href="fixes.html">Fix Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="computes.html">Compute Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="pairs.html">Pair Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="bonds.html">Bond Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="angles.html">Angle Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="dihedrals.html">Dihedral Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="impropers.html">Improper Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="dumps.html">Dump Styles</a></li>
<li class="toctree-l1"><a class="reference internal" href="fix_modify_atc_commands.html">fix_modify AtC commands</a></li>
<li class="toctree-l1"><a class="reference internal" href="Bibliography.html">Bibliography</a></li>
</ul>
</div>
</div>
</nav>
<section data-toggle="wy-nav-shift" class="wy-nav-content-wrap"><nav class="wy-nav-top" aria-label="Mobile navigation menu" >
<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 style-external-links">
<div role="navigation" aria-label="Page navigation">
<ul class="wy-breadcrumbs">
<li><a href="Manual.html" class="icon icon-home" aria-label="Home"></a></li>
<li class="breadcrumb-item"><a href="Howto.html"><span class="section-number">8. </span>Howto discussions</a></li>
<li class="breadcrumb-item active"><span class="section-number">8.5.6. </span>Adiabatic core/shell model</li>
<li class="wy-breadcrumbs-aside">
<a href="https://www.lammps.org"><img src="_static/lammps-logo.png" width="64" height="16" alt="LAMMPS Homepage"></a> | <a href="Commands_all.html">Commands</a>
</li>
</ul><div class="rst-breadcrumbs-buttons" role="navigation" aria-label="Sequential page navigation">
<a href="Howto_polarizable.html" class="btn btn-neutral float-left" title="8.5.5. Polarizable models" accesskey="p"><span class="fa fa-arrow-circle-left" aria-hidden="true"></span> Previous</a>
<a href="Howto_drude.html" class="btn btn-neutral float-right" title="8.5.7. Drude induced dipoles" accesskey="n">Next <span class="fa fa-arrow-circle-right" aria-hidden="true"></span></a>
</div>
<hr/>
</div>
<div role="main" class="document" itemscope="itemscope" itemtype="http://schema.org/Article">
<div itemprop="articleBody">
<p><span class="math notranslate nohighlight">\(\renewcommand{\AA}{\text{Å}}\)</span></p>
<section id="adiabatic-core-shell-model">
<h1><span class="section-number">8.5.6. </span>Adiabatic core/shell model<a class="headerlink" href="#adiabatic-core-shell-model" title="Link to this heading"></a></h1>
<p>The adiabatic core-shell model by <a class="reference internal" href="#mitchellfincham"><span class="std std-ref">Mitchell and Fincham</span></a> is a simple method for adding polarizability
to a system. In order to mimic the electron shell of an ion, a
satellite particle is attached to it. This way the ions are split into
a core and a shell where the latter is meant to react to the
electrostatic environment inducing polarizability. See the <a class="reference internal" href="Howto_polarizable.html"><span class="doc">Howto polarizable</span></a> page for a discussion of all
the polarizable models available in LAMMPS.</p>
<p>Technically, shells are attached to the cores by a spring force f =
k*r where k is a parameterized spring constant and r is the distance
between the core and the shell. The charges of the core and the shell
add up to the ion charge, thus q(ion) = q(core) + q(shell). This
setup introduces the ion polarizability (alpha) given by
alpha = q(shell)^2 / k. In a
similar fashion the mass of the ion is distributed on the core and the
shell with the core having the larger mass.</p>
<p>To run this model in LAMMPS, <a class="reference internal" href="atom_style.html"><span class="doc">atom_style</span></a> <em>full</em> can
be used since atom charge and bonds are needed. Each kind of
core/shell pair requires two atom types and a bond type. The core and
shell of a core/shell pair should be bonded to each other with a
harmonic bond that provides the spring force. For example, a data file
for NaCl, as found in examples/coreshell, has this format:</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>432 atoms # core and shell atoms
216 bonds # number of core/shell springs
4 atom types # 2 cores and 2 shells for Na and Cl
2 bond types
0.0 24.09597 xlo xhi
0.0 24.09597 ylo yhi
0.0 24.09597 zlo zhi
Masses # core/shell mass ratio = 0.1
1 20.690784 # Na core
2 31.90500 # Cl core
3 2.298976 # Na shell
4 3.54500 # Cl shell
Atoms
1 1 2 1.5005 0.00000000 0.00000000 0.00000000 # core of core/shell pair 1
2 1 4 -2.5005 0.00000000 0.00000000 0.00000000 # shell of core/shell pair 1
3 2 1 1.5056 4.01599500 4.01599500 4.01599500 # core of core/shell pair 2
4 2 3 -0.5056 4.01599500 4.01599500 4.01599500 # shell of core/shell pair 2
(...)
Bonds # Bond topology for spring forces
1 2 1 2 # spring for core/shell pair 1
2 2 3 4 # spring for core/shell pair 2
(...)
</pre></div>
</div>
<p>Non-Coulombic (e.g. Lennard-Jones) pairwise interactions are only
defined between the shells. Coulombic interactions are defined
between all cores and shells. If desired, additional bonds can be
specified between cores.</p>
<p>The <a class="reference internal" href="special_bonds.html"><span class="doc">special_bonds</span></a> command should be used to
turn-off the Coulombic interaction within core/shell pairs, since that
interaction is set by the bond spring. This is done using the
<a class="reference internal" href="special_bonds.html"><span class="doc">special_bonds</span></a> command with a 1-2 weight = 0.0,
which is the default value. It needs to be considered whether one has
to adjust the <a class="reference internal" href="special_bonds.html"><span class="doc">special_bonds</span></a> weighting according
to the molecular topology since the interactions of the shells are
bypassed over an extra bond.</p>
<p>Note that this core/shell implementation does not require all ions to
be polarized. One can mix core/shell pairs and ions without a
satellite particle if desired.</p>
<p>Since the core/shell model permits distances of r = 0.0 between the
core and shell, a pair style with a “cs” suffix needs to be used to
implement a valid long-range Coulombic correction. Several such pair
styles are provided in the CORESHELL package. See <a class="reference internal" href="pair_cs.html"><span class="doc">this page</span></a> for details. All of the core/shell enabled pair
styles require the use of a long-range Coulombic solver, as specified
by the <a class="reference internal" href="kspace_style.html"><span class="doc">kspace_style</span></a> command. Either the PPPM or
Ewald solvers can be used.</p>
<p>For the NaCL example problem, these pair style and bond style settings
are used:</p>
<div class="highlight-LAMMPS notranslate"><div class="highlight"><pre><span></span><span class="k">pair_style</span><span class="w"> </span><span class="n">born</span><span class="o">/</span><span class="n">coul</span><span class="o">/</span><span class="n">long</span><span class="o">/</span><span class="n">cs</span><span class="w"> </span><span class="m">20.0</span><span class="w"> </span><span class="m">20.0</span>
<span class="k">pair_coeff</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="o">*</span><span class="w"> </span><span class="m">0.0</span><span class="w"> </span><span class="m">1.000</span><span class="w"> </span><span class="m">0.00</span><span class="w"> </span><span class="m">0.00</span><span class="w"> </span><span class="m">0.00</span>
<span class="k">pair_coeff</span><span class="w"> </span><span class="m">3</span><span class="w"> </span><span class="m">3</span><span class="w"> </span><span class="m">487.0</span><span class="w"> </span><span class="m">0.23768</span><span class="w"> </span><span class="m">0.00</span><span class="w"> </span><span class="m">1.05</span><span class="w"> </span><span class="m">0.50</span><span class="w"> </span><span class="c">#Na-Na</span>
<span class="k">pair_coeff</span><span class="w"> </span><span class="m">3</span><span class="w"> </span><span class="m">4</span><span class="w"> </span><span class="m">145134.0</span><span class="w"> </span><span class="m">0.23768</span><span class="w"> </span><span class="m">0.00</span><span class="w"> </span><span class="m">6.99</span><span class="w"> </span><span class="m">8.70</span><span class="w"> </span><span class="c">#Na-Cl</span>
<span class="k">pair_coeff</span><span class="w"> </span><span class="m">4</span><span class="w"> </span><span class="m">4</span><span class="w"> </span><span class="m">405774.0</span><span class="w"> </span><span class="m">0.23768</span><span class="w"> </span><span class="m">0.00</span><span class="w"> </span><span class="m">72.40</span><span class="w"> </span><span class="m">145.40</span><span class="w"> </span><span class="c">#Cl-Cl</span>
<span class="k">bond_style</span><span class="w"> </span><span class="n">harmonic</span>
<span class="k">bond_coeff</span><span class="w"> </span><span class="m">1</span><span class="w"> </span><span class="m">63.014</span><span class="w"> </span><span class="m">0.0</span>
<span class="k">bond_coeff</span><span class="w"> </span><span class="m">2</span><span class="w"> </span><span class="m">25.724</span><span class="w"> </span><span class="m">0.0</span>
</pre></div>
</div>
<p>When running dynamics with the adiabatic core/shell model, the
following issues should be considered. The relative motion of
the core and shell particles corresponds to the polarization,
hereby an instantaneous relaxation of the shells is approximated
and a fast core/shell spring frequency ensures a nearly constant
internal kinetic energy during the simulation.
Thermostats can alter this polarization behavior, by scaling the
internal kinetic energy, meaning the shell will not react freely to
its electrostatic environment.
Therefore it is typically desirable to decouple the relative motion of
the core/shell pair, which is an imaginary degree of freedom, from the
real physical system. To do that, the <a class="reference internal" href="compute_temp_cs.html"><span class="doc">compute temp/cs</span></a> command can be used, in conjunction with
any of the thermostat fixes, such as <a class="reference internal" href="fix_nh.html"><span class="doc">fix nvt</span></a> or <a class="reference internal" href="fix_langevin.html"><span class="doc">fix langevin</span></a>. This compute uses the center-of-mass velocity
of the core/shell pairs to calculate a temperature, and ensures that
velocity is what is rescaled for thermostatting purposes. This
compute also works for a system with both core/shell pairs and
non-polarized ions (ions without an attached satellite particle). The
<a class="reference internal" href="compute_temp_cs.html"><span class="doc">compute temp/cs</span></a> command requires input of two
groups, one for the core atoms, another for the shell atoms.
Non-polarized ions which might also be included in the treated system
should not be included into either of these groups, they are taken
into account by the <em>group-ID</em> (second argument) of the compute. The
groups can be defined using the <a class="reference internal" href="group.html"><span class="doc">group *type*</span></a> command.
Note that to perform thermostatting using this definition of
temperature, the <a class="reference internal" href="fix_modify.html"><span class="doc">fix modify temp</span></a> command should be
used to assign the compute to the thermostat fix. Likewise the
<a class="reference internal" href="thermo_modify.html"><span class="doc">thermo_modify temp</span></a> command can be used to make
this temperature be output for the overall system.</p>
<p>For the NaCl example, this can be done as follows:</p>
<div class="highlight-LAMMPS notranslate"><div class="highlight"><pre><span></span><span class="k">group </span><span class="nv nv-Identifier">cores</span><span class="w"> </span><span class="n">type</span><span class="w"> </span><span class="m">1</span><span class="w"> </span><span class="m">2</span>
<span class="k">group </span><span class="nv nv-Identifier">shells</span><span class="w"> </span><span class="n">type</span><span class="w"> </span><span class="m">3</span><span class="w"> </span><span class="m">4</span>
<span class="k">compute </span><span class="nv nv-Identifier">CSequ</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">temp</span><span class="o">/</span><span class="n">cs</span><span class="w"> </span><span class="n">cores</span><span class="w"> </span><span class="n">shells</span>
<span class="k">fix </span><span class="nv nv-Identifier">thermoberendsen</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">temp</span><span class="o">/</span><span class="n">berendsen</span><span class="w"> </span><span class="m">1427</span><span class="w"> </span><span class="m">1427</span><span class="w"> </span><span class="m">0.4</span><span class="w"> </span><span class="c"># thermostat for the true physical system</span>
<span class="k">fix </span><span class="nv nv-Identifier">thermostatequ</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">nve</span><span class="w"> </span><span class="c"># integrator as needed for the berendsen thermostat</span>
<span class="k">fix_modify </span><span class="nv nv-Identifier">thermoberendsen</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">CSequ</span>
<span class="k">thermo_modify</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">CSequ</span><span class="w"> </span><span class="c"># output of center-of-mass derived temperature</span>
</pre></div>
</div>
<p>The pressure for the core/shell system is computed via the regular
LAMMPS convention by <a class="reference internal" href="#mitchellfincham2"><span class="std std-ref">treating the cores and shells as individual particles</span></a>. For the thermo output of the pressure
as well as for the application of a barostat, it is necessary to
use an additional <a class="reference internal" href="compute_pressure.html"><span class="doc">pressure</span></a> compute based on
the default <a class="reference internal" href="compute_temp.html"><span class="doc">temperature</span></a> and specifying it as a
second argument in <a class="reference internal" href="fix_modify.html"><span class="doc">fix modify</span></a> and
<a class="reference internal" href="thermo_modify.html"><span class="doc">thermo_modify</span></a> resulting in:</p>
<div class="highlight-LAMMPS notranslate"><div class="highlight"><pre><span></span><span class="nv">(...)</span>
<span class="k">compute </span><span class="nv nv-Identifier">CSequ</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">temp</span><span class="o">/</span><span class="n">cs</span><span class="w"> </span><span class="n">cores</span><span class="w"> </span><span class="n">shells</span>
<span class="k">compute </span><span class="nv nv-Identifier">thermo_press_lmp</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">pressure</span><span class="w"> </span><span class="n">thermo_temp</span><span class="w"> </span><span class="c"># pressure for individual particles</span>
<span class="k">thermo_modify</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">CSequ</span><span class="w"> </span><span class="n">press</span><span class="w"> </span><span class="n">thermo_press_lmp</span><span class="w"> </span><span class="c"># modify thermo to regular pressure</span>
<span class="k">fix </span><span class="nv nv-Identifier">press_bar</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">npt</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="m">300</span><span class="w"> </span><span class="m">300</span><span class="w"> </span><span class="m">0.04</span><span class="w"> </span><span class="n">iso</span><span class="w"> </span><span class="m">0</span><span class="w"> </span><span class="m">0</span><span class="w"> </span><span class="m">0.4</span>
<span class="k">fix_modify </span><span class="nv nv-Identifier">press_bar</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">CSequ</span><span class="w"> </span><span class="n">press</span><span class="w"> </span><span class="n">thermo_press_lmp</span><span class="w"> </span><span class="c"># pressure modification for correct kinetic scalar</span>
</pre></div>
</div>
<p>If <a class="reference internal" href="compute_temp_cs.html"><span class="doc">compute temp/cs</span></a> is used, the decoupled
relative motion of the core and the shell should in theory be
stable. However numerical fluctuation can introduce a small
momentum to the system, which is noticeable over long trajectories.
Therefore it is recommendable to use the <a class="reference internal" href="fix_momentum.html"><span class="doc">fix momentum</span></a> command in combination with <a class="reference internal" href="compute_temp_cs.html"><span class="doc">compute temp/cs</span></a> when equilibrating the system to
prevent any drift.</p>
<p>When initializing the velocities of a system with core/shell pairs, it
is also desirable to not introduce energy into the relative motion of
the core/shell particles, but only assign a center-of-mass velocity to
the pairs. This can be done by using the <em>bias</em> keyword of the
<a class="reference internal" href="velocity.html"><span class="doc">velocity create</span></a> command and assigning the <a class="reference internal" href="compute_temp_cs.html"><span class="doc">compute temp/cs</span></a> command to the <em>temp</em> keyword of the
<a class="reference internal" href="velocity.html"><span class="doc">velocity</span></a> command, e.g.</p>
<div class="highlight-LAMMPS notranslate"><div class="highlight"><pre><span></span><span class="k">velocity </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">create</span><span class="w"> </span><span class="m">1427</span><span class="w"> </span><span class="m">134</span><span class="w"> </span><span class="n">bias</span><span class="w"> </span><span class="n">yes</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">CSequ</span>
<span class="k">velocity </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">scale</span><span class="w"> </span><span class="m">1427</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">CSequ</span>
</pre></div>
</div>
<p>To maintain the correct polarizability of the core/shell pairs, the
kinetic energy of the internal motion shall remain nearly constant.
Therefore the choice of spring force and mass ratio need to ensure
much faster relative motion of the two atoms within the core/shell pair
than their center-of-mass velocity. This allows the shells to
effectively react instantaneously to the electrostatic environment and
limits energy transfer to or from the core/shell oscillators.
This fast movement also dictates the timestep that can be used.</p>
<p>The primary literature of the adiabatic core/shell model suggests that
the fast relative motion of the core/shell pairs only allows negligible
energy transfer to the environment.
The mentioned energy transfer will typically lead to a small drift
in total energy over time. This internal energy can be monitored
using the <a class="reference internal" href="compute_chunk_atom.html"><span class="doc">compute chunk/atom</span></a> and <a class="reference internal" href="compute_temp_chunk.html"><span class="doc">compute temp/chunk</span></a> commands. The internal kinetic
energies of each core/shell pair can then be summed using the sum()
special function of the <a class="reference internal" href="variable.html"><span class="doc">variable</span></a> command. Or they can
be time/averaged and output using the <a class="reference internal" href="fix_ave_time.html"><span class="doc">fix ave/time</span></a>
command. To use these commands, each core/shell pair must be defined
as a “chunk”. If each core/shell pair is defined as its own molecule,
the molecule ID can be used to define the chunks. If cores are bonded
to each other to form larger molecules, the chunks can be identified
by the <a class="reference internal" href="fix_property_atom.html"><span class="doc">fix property/atom</span></a> via assigning a
core/shell ID to each atom using a special field in the data file read
by the <a class="reference internal" href="read_data.html"><span class="doc">read_data</span></a> command. This field can then be
accessed by the <a class="reference internal" href="compute_property_atom.html"><span class="doc">compute property/atom</span></a>
command, to use as input to the <a class="reference internal" href="compute_chunk_atom.html"><span class="doc">compute chunk/atom</span></a> command to define the core/shell
pairs as chunks.</p>
<p>For example if core/shell pairs are the only molecules:</p>
<div class="highlight-LAMMPS notranslate"><div class="highlight"><pre><span></span><span class="k">read_data</span><span class="w"> </span><span class="n">NaCl_CS_x0.1_prop.data</span>
<span class="k">compute </span><span class="nv nv-Identifier">prop</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">property</span><span class="o">/</span><span class="n">atom</span><span class="w"> </span><span class="n">molecule</span>
<span class="k">compute </span><span class="nv nv-Identifier">cs_chunk</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">chunk</span><span class="o">/</span><span class="n">atom</span><span class="w"> </span><span class="n">c_prop</span>
<span class="k">compute </span><span class="nv nv-Identifier">cstherm</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">temp</span><span class="o">/</span><span class="n">chunk</span><span class="w"> </span><span class="n">cs_chunk</span><span class="w"> </span><span class="n">temp</span><span class="w"> </span><span class="n">internal</span><span class="w"> </span><span class="n">com</span><span class="w"> </span><span class="n">yes</span><span class="w"> </span><span class="n">cdof</span><span class="w"> </span><span class="m">3.0</span><span class="w"> </span><span class="c"># note the chosen degrees of freedom for the core/shell pairs</span>
<span class="k">fix </span><span class="nv nv-Identifier">ave_chunk</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">ave</span><span class="o">/</span><span class="n">time</span><span class="w"> </span><span class="m">10</span><span class="w"> </span><span class="m">1</span><span class="w"> </span><span class="m">10</span><span class="w"> </span><span class="n">c_cstherm</span><span class="w"> </span><span class="n">file</span><span class="w"> </span><span class="n">chunk.dump</span><span class="w"> </span><span class="n">mode</span><span class="w"> </span><span class="n">vector</span>
</pre></div>
</div>
<p>For example if core/shell pairs and other molecules are present:</p>
<div class="highlight-LAMMPS notranslate"><div class="highlight"><pre><span></span><span class="k">fix </span><span class="nv nv-Identifier">csinfo</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">property</span><span class="o">/</span><span class="n">atom</span><span class="w"> </span><span class="n">i_CSID</span><span class="w"> </span><span class="c"># property/atom command</span>
<span class="k">read_data</span><span class="w"> </span><span class="n">NaCl_CS_x0.1_prop.data</span><span class="w"> </span><span class="k">fix </span><span class="nv nv-Identifier">csinfo</span><span class="w"> </span><span class="nv nv-Identifier">NULL</span><span class="w"> </span><span class="n">CS</span><span class="o">-</span><span class="n">Info</span><span class="w"> </span><span class="c"># atom property added in the data-file</span>
<span class="k">compute </span><span class="nv nv-Identifier">prop</span><span class="w"> </span><span class="nv nv-Identifier">all</span><span class="w"> </span><span class="n">property</span><span class="o">/</span><span class="n">atom</span><span class="w"> </span><span class="n">i_CSID</span>
<span class="nv">(...)</span>
</pre></div>
</div>
<p>The additional section in the date file would be formatted like this:</p>
<div class="highlight-none notranslate"><div class="highlight"><pre><span></span>CS-Info # header of additional section
1 1 # column 1 = atom ID, column 2 = core/shell ID
2 1
3 2
4 2
5 3
6 3
7 4
8 4
(...)
</pre></div>
</div>
<hr class="docutils" />
<p id="mitchellfincham"><strong>(Mitchell and Fincham)</strong> Mitchell, Fincham, J Phys Condensed Matter,
5, 1031-1038 (1993).</p>
<p id="mitchellfincham2"><strong>(Fincham)</strong> Fincham, Mackrodt and Mitchell, J Phys Condensed Matter,
6, 393-404 (1994).</p>
</section>
</div>
</div>
<footer><div class="rst-footer-buttons" role="navigation" aria-label="Footer">
<a href="Howto_polarizable.html" class="btn btn-neutral float-left" title="8.5.5. Polarizable models" accesskey="p" rel="prev"><span class="fa fa-arrow-circle-left" aria-hidden="true"></span> Previous</a>
<a href="Howto_drude.html" class="btn btn-neutral float-right" title="8.5.7. Drude induced dipoles" accesskey="n" rel="next">Next <span class="fa fa-arrow-circle-right" aria-hidden="true"></span></a>
</div>
<hr/>
<div role="contentinfo">
<p>&#169; Copyright 2003-2025 Sandia Corporation.</p>
</div>
Built with <a href="https://www.sphinx-doc.org/">Sphinx</a> using a
<a href="https://github.com/readthedocs/sphinx_rtd_theme">theme</a>
provided by <a href="https://readthedocs.org">Read the Docs</a>.
</footer>
</div>
</div>
</section>
</div>
<script>
jQuery(function () {
SphinxRtdTheme.Navigation.enable(false);
});
</script>
</body>
</html>