mirror of https://github.com/lammps/lammps.git
89 lines
3.4 KiB
Plaintext
89 lines
3.4 KiB
Plaintext
# file "tip3p2004.lt"
|
|
#
|
|
# H1 H2
|
|
# \ /
|
|
# O
|
|
#
|
|
# I think this is the TIP3P water described in the paper by
|
|
# Daniel J. Price and Charles L. Brooks III
|
|
# J. Chem. Phys., 121(20): 10096 (2004)
|
|
# Specifically I think it refers to the "Model B" version of long-range TIP3P
|
|
# described in the 3rd-to-last column of "Table I", on p.10099.
|
|
|
|
TIP3P2004 {
|
|
|
|
write_once("In Init") {
|
|
# -- Default styles (for solo "TIP3P2004" water) --
|
|
units real
|
|
atom_style full
|
|
|
|
pair_style hybrid lj/charmm/coul/long 10.0 10.5 10.5
|
|
bond_style hybrid harmonic
|
|
angle_style hybrid harmonic
|
|
kspace_style pppm 0.0001
|
|
pair_modify mix arithmetic
|
|
}
|
|
|
|
write("Data Atoms") {
|
|
$atom:O $mol:. @atom:O -0.830 0.0000000 0.00000 0.000000
|
|
$atom:H1 $mol:. @atom:H 0.415 0.756950327 0.00000 0.5858822766
|
|
$atom:H2 $mol:. @atom:H 0.415 -0.756950327 0.00000 0.5858822766
|
|
}
|
|
|
|
write_once("Data Masses") {
|
|
@atom:O 15.9994
|
|
@atom:H 1.008
|
|
}
|
|
|
|
write("Data Bonds") {
|
|
$bond:OH1 @bond:OH $atom:O $atom:H1
|
|
$bond:OH2 @bond:OH $atom:O $atom:H2
|
|
}
|
|
|
|
write("Data Angles") {
|
|
$angle:HOH @angle:HOH $atom:H1 $atom:O $atom:H2
|
|
}
|
|
|
|
write_once("In Settings") {
|
|
bond_coeff @bond:OH harmonic 450.0 0.9572
|
|
angle_coeff @angle:HOH harmonic 55.0 104.52
|
|
|
|
#########################################################################
|
|
#### There are two choices for for the O-O interactions
|
|
#########################################################################
|
|
#### O-O nonbonded interactions
|
|
# For the 1983 Jorgensen version of TIP3P use:
|
|
# pair_coeff @atom:O @atom:O lj/charmm/coul/charmm 0.1521 3.1507
|
|
# For the 2004 Price & Brooks version of TIP3P use:
|
|
pair_coeff @atom:O @atom:O lj/charmm/coul/charmm 0.102 3.188
|
|
#########################################################################
|
|
#### There are three choices for for the O-H and H-H interactions
|
|
#########################################################################
|
|
#### 1) CHARMM uses an arithmetic mixing-rule for the O-H sigma parameter
|
|
pair_coeff @atom:H @atom:H lj/charmm/coul/charmm 0.0460 0.4000
|
|
pair_coeff @atom:O @atom:H lj/charmm/coul/charmm 0.0836 1.7753 #arithmetic
|
|
#########################################################################
|
|
#### 2) OPLS-AA uses geometric a mixing-fule for the O-H sigma parameter,
|
|
#### If you want to use this, uncomment the following two lines:
|
|
# pair_coeff @atom:H @atom:H lj/charmm/coul/charmm 0.0460 0.4000
|
|
# pair_coeff @atom:O @atom:H lj/charmm/coul/charmm 0.0836 1.1226 #geometric
|
|
#########################################################################
|
|
#### 3) The original Jorgensen 1983 parameterization has no OH or HH
|
|
# lennard-jones interactions. For this behavior, uncomment these lines:
|
|
# pair_coeff @atom:H @atom:H lj/charmm/coul/charmm 0.00 0.4000
|
|
# pair_coeff @atom:O @atom:H lj/charmm/coul/charmm 0.00 1.7753
|
|
#########################################################################
|
|
|
|
# Define a group for the tip3p water molecules:
|
|
group tip3p type @atom:O @atom:H
|
|
|
|
# Optional: Constrain the angles and distances.
|
|
# (Most implementations use this, but it is optional.)
|
|
fix fSHAKE tip3p shake 0.0001 10 100 b @bond:OH a @angle:HOH
|
|
# (Remember to "unfix" fSHAKE during minimization.)
|
|
|
|
}
|
|
|
|
} # "TIP3P2004" water molecule type
|
|
|