"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

pair_style oxdna/excv command :h3
pair_style oxdna/stk command :h3
pair_style oxdna/hbond command :h3
pair_style oxdna/xstk command :h3
pair_style oxdna/coaxstk command :h3

[Syntax:]

pair_style style1 :pre

pair_coeff * * style2 args :pre

style1 = {hybrid/overlay oxdna/excv oxdna/stk oxdna/hbond oxdna/xstk oxdna/coaxstk} :ul

style2 = {oxdna/excv} or {oxdna/stk} or {oxdna/hbond} or {oxdna/xstk} or {oxdna/coaxstk}
args = list of arguments for these particular styles :ul

  {oxdna/stk} args = seq T 6.0 0.4 0.9 0.32 0.6 1.3 0 0.8 0.9 0 0.95 0.9 0 0.95 2.0 0.65 2.0 0.65
    seq = seqav (for average sequence stacking strength) or seqdep (for sequence-dependent stacking strength)
    T = temperature (oxDNA units, 0.1 = 300 K) 
  {oxdna/hbond} args = seq eps 8.0 0.4 0.75 0.34 0.7 1.5 0 0.7 1.5 0 0.7 1.5 0 0.7 0.46 3.141592653589793 0.7 4.0 1.5707963267948966 0.45 4.0 1.5707963267948966 0.45
    seq = seqav (for average sequence base-pairing strength) or seqdep (for sequence-dependent base-pairing strength)
    eps = 1.077 (between base pairs A-T and C-G) or 0 (all other pairs) :pre 

[Examples:]

pair_style hybrid/overlay oxdna/excv oxdna/stk oxdna/hbond oxdna/xstk oxdna/coaxstk
pair_coeff * * oxdna/excv    2.0 0.7 0.675 2.0 0.515 0.5 2.0 0.33 0.32
pair_coeff * * oxdna/stk     seqdep 0.1 6.0 0.4 0.9 0.32 0.6 1.3 0 0.8 0.9 0 0.95 0.9 0 0.95 2.0 0.65 2.0 0.65
pair_coeff * * oxdna/hbond   seqdep 0.0 8.0 0.4 0.75 0.34 0.7 1.5 0 0.7 1.5 0 0.7 1.5 0 0.7 0.46 3.141592653589793 0.7 4.0 1.5707963267948966 0.45 4.0 1.5707963267948966 0.45
pair_coeff 1 4 oxdna/hbond   seqdep 1.077 8.0 0.4 0.75 0.34 0.7 1.5 0 0.7 1.5 0 0.7 1.5 0 0.7 0.46 3.141592653589793 0.7 4.0 1.5707963267948966 0.45 4.0 1.5707963267948966 0.45
pair_coeff 2 3 oxdna/hbond   seqdep 1.077 8.0 0.4 0.75 0.34 0.7 1.5 0 0.7 1.5 0 0.7 1.5 0 0.7 0.46 3.141592653589793 0.7 4.0 1.5707963267948966 0.45 4.0 1.5707963267948966 0.45
pair_coeff * * oxdna/xstk    47.5 0.575 0.675 0.495 0.655 2.25 0.791592653589793 0.58 1.7 1.0 0.68 1.7 1.0 0.68 1.5 0 0.65 1.7 0.875 0.68 1.7 0.875 0.68
pair_coeff * * oxdna/coaxstk 46.0 0.4 0.6 0.22 0.58 2.0 2.541592653589793 0.65 1.3 0 0.8 0.9 0 0.95 0.9 0 0.95 2.0 -0.65 2.0 -0.65 :pre

[Description:]

The {oxdna} pair styles compute the pairwise-additive parts of the oxDNA force field
for coarse-grained modelling of DNA. The effective interaction between the nucleotides consists of potentials for the
excluded volume interaction {oxdna/excv}, the stacking {oxdna/stk}, cross-stacking {oxdna/xstk}
and coaxial stacking interaction {oxdna/coaxstk} as well
as the hydrogen-bonding interaction {oxdna/hbond} between complementary pairs of nucleotides on
opposite strands. Average sequence or sequence-dependent stacking and base-pairing strengths
are supported "(Sulc)"_#Sulc1.

The exact functional form of the pair styles is rather complex.
The individual potentials consist of products of modulation factors,
which themselves are constructed from a number of more basic potentials
(Morse, Lennard-Jones, harmonic angle and distance) as well as quadratic smoothing and modulation terms.
We refer to "(Ouldridge-DPhil)"_#Ouldridge-DPhil1 and "(Ouldridge)"_#Ouldridge1
for a detailed description of the oxDNA force field.

NOTE: These pair styles have to be used together with the related oxDNA bond style
{oxdna/fene} for the connectivity of the phosphate backbone (see also documentation of
"bond_style oxdna/fene"_bond_oxdna.html). Most of the coefficients
in the above example have to be kept fixed and cannot be changed without reparametrizing the entire model.
Exceptions are the first and second coefficient after {oxdna/stk} (seq=seqdep and T=0.1 in the above example)
and the first coefficient after {oxdna/hbond} (seq=seqdep in the above example).
When using a Langevin thermostat, e.g. through "fix langevin"_fix_langevin.html
or "fix nve/dotc/langevin"_fix_nve_dotc_langevin.html
the temperature coefficients have to be matched to the one used in the fix.

Example input and data files for DNA duplexes can be found in examples/USER/cgdna/examples/oxDNA/ and /oxDNA2/.
A simple python setup tool which creates single straight or helical DNA strands,
DNA duplexes or arrays of DNA duplexes can be found in examples/USER/cgdna/util/.
A technical report with more information on the model, the structure of the input file,
the setup tool and the performance of the LAMMPS-implementation of oxDNA
can be found "here"_PDF/USER-CGDNA-overview.pdf.

:line

[Restrictions:]

These pair styles can only be used if LAMMPS was built with the
USER-CGDNA package and the MOLECULE and ASPHERE package.  See the "Making
LAMMPS"_Section_start.html#start_3 section for more info on packages.

[Related commands:]

"bond_style oxdna/fene"_bond_oxdna.html, "fix nve/dotc/langevin"_fix_nve_dotc_langevin.html, "pair_coeff"_pair_coeff.html,
"bond_style oxdna2/fene"_bond_oxdna.html, "pair_style oxdna2/excv"_pair_oxdna2.html

[Default:] none

:line

:link(Sulc1)
[(Sulc)] P. Sulc, F. Romano, T.E. Ouldridge, L. Rovigatti, J.P.K. Doye, A.A. Louis, J. Chem. Phys. 137, 135101 (2012).

:link(Ouldridge-DPhil1)
[(Ouldrigde-DPhil)] T.E. Ouldridge, Coarse-grained modelling of DNA and DNA self-assembly, DPhil. University of Oxford (2011).

:link(Ouldridge1)
[(Ouldridge)] T.E. Ouldridge, A.A. Louis, J.P.K. Doye, J. Chem. Phys. 134, 085101 (2011).
