Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F91303786
pair_lj_class2_coul_long.cpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sat, Nov 9, 20:06
Size
13 KB
Mime Type
text/x-c
Expires
Mon, Nov 11, 20:06 (2 d)
Engine
blob
Format
Raw Data
Handle
22239590
Attached To
rLAMMPS lammps
pair_lj_class2_coul_long.cpp
View Options
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#include "math.h"
#include "stdio.h"
#include "stdlib.h"
#include "string.h"
#include "pair_lj_class2_coul_long.h"
#include "atom.h"
#include "comm.h"
#include "force.h"
#include "kspace.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "math_const.h"
#include "memory.h"
#include "error.h"
using
namespace
LAMMPS_NS
;
using
namespace
MathConst
;
#define EWALD_F 1.12837917
#define EWALD_P 0.3275911
#define A1 0.254829592
#define A2 -0.284496736
#define A3 1.421413741
#define A4 -1.453152027
#define A5 1.061405429
/* ---------------------------------------------------------------------- */
PairLJClass2CoulLong
::
PairLJClass2CoulLong
(
LAMMPS
*
lmp
)
:
Pair
(
lmp
)
{}
/* ---------------------------------------------------------------------- */
PairLJClass2CoulLong
::~
PairLJClass2CoulLong
()
{
if
(
allocated
)
{
memory
->
destroy
(
setflag
);
memory
->
destroy
(
cutsq
);
memory
->
destroy
(
cut_lj
);
memory
->
destroy
(
cut_ljsq
);
memory
->
destroy
(
epsilon
);
memory
->
destroy
(
sigma
);
memory
->
destroy
(
lj1
);
memory
->
destroy
(
lj2
);
memory
->
destroy
(
lj3
);
memory
->
destroy
(
lj4
);
memory
->
destroy
(
offset
);
}
}
/* ---------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
compute
(
int
eflag
,
int
vflag
)
{
int
i
,
j
,
ii
,
jj
,
inum
,
jnum
,
itype
,
jtype
;
double
qtmp
,
xtmp
,
ytmp
,
ztmp
,
delx
,
dely
,
delz
,
evdwl
,
ecoul
,
fpair
;
double
rsq
,
r
,
rinv
,
r2inv
,
r3inv
,
r6inv
,
forcecoul
,
forcelj
;
double
grij
,
expm2
,
prefactor
,
t
,
erfc
;
double
factor_coul
,
factor_lj
;
int
*
ilist
,
*
jlist
,
*
numneigh
,
**
firstneigh
;
evdwl
=
ecoul
=
0.0
;
if
(
eflag
||
vflag
)
ev_setup
(
eflag
,
vflag
);
else
evflag
=
vflag_fdotr
=
0
;
double
**
x
=
atom
->
x
;
double
**
f
=
atom
->
f
;
double
*
q
=
atom
->
q
;
int
*
type
=
atom
->
type
;
int
nlocal
=
atom
->
nlocal
;
double
*
special_coul
=
force
->
special_coul
;
double
*
special_lj
=
force
->
special_lj
;
int
newton_pair
=
force
->
newton_pair
;
double
qqrd2e
=
force
->
qqrd2e
;
inum
=
list
->
inum
;
ilist
=
list
->
ilist
;
numneigh
=
list
->
numneigh
;
firstneigh
=
list
->
firstneigh
;
// loop over neighbors of my atoms
for
(
ii
=
0
;
ii
<
inum
;
ii
++
)
{
i
=
ilist
[
ii
];
qtmp
=
q
[
i
];
xtmp
=
x
[
i
][
0
];
ytmp
=
x
[
i
][
1
];
ztmp
=
x
[
i
][
2
];
itype
=
type
[
i
];
jlist
=
firstneigh
[
i
];
jnum
=
numneigh
[
i
];
for
(
jj
=
0
;
jj
<
jnum
;
jj
++
)
{
j
=
jlist
[
jj
];
factor_lj
=
special_lj
[
sbmask
(
j
)];
factor_coul
=
special_coul
[
sbmask
(
j
)];
j
&=
NEIGHMASK
;
delx
=
xtmp
-
x
[
j
][
0
];
dely
=
ytmp
-
x
[
j
][
1
];
delz
=
ztmp
-
x
[
j
][
2
];
rsq
=
delx
*
delx
+
dely
*
dely
+
delz
*
delz
;
jtype
=
type
[
j
];
if
(
rsq
<
cutsq
[
itype
][
jtype
])
{
r2inv
=
1.0
/
rsq
;
if
(
rsq
<
cut_coulsq
)
{
r
=
sqrt
(
rsq
);
grij
=
g_ewald
*
r
;
expm2
=
exp
(
-
grij
*
grij
);
t
=
1.0
/
(
1.0
+
EWALD_P
*
grij
);
erfc
=
t
*
(
A1
+
t
*
(
A2
+
t
*
(
A3
+
t
*
(
A4
+
t
*
A5
))))
*
expm2
;
prefactor
=
qqrd2e
*
qtmp
*
q
[
j
]
/
r
;
forcecoul
=
prefactor
*
(
erfc
+
EWALD_F
*
grij
*
expm2
);
if
(
factor_coul
<
1.0
)
forcecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
}
else
forcecoul
=
0.0
;
if
(
rsq
<
cut_ljsq
[
itype
][
jtype
])
{
rinv
=
sqrt
(
r2inv
);
r3inv
=
r2inv
*
rinv
;
r6inv
=
r3inv
*
r3inv
;
forcelj
=
r6inv
*
(
lj1
[
itype
][
jtype
]
*
r3inv
-
lj2
[
itype
][
jtype
]);
}
else
forcelj
=
0.0
;
fpair
=
(
forcecoul
+
factor_lj
*
forcelj
)
*
r2inv
;
f
[
i
][
0
]
+=
delx
*
fpair
;
f
[
i
][
1
]
+=
dely
*
fpair
;
f
[
i
][
2
]
+=
delz
*
fpair
;
if
(
newton_pair
||
j
<
nlocal
)
{
f
[
j
][
0
]
-=
delx
*
fpair
;
f
[
j
][
1
]
-=
dely
*
fpair
;
f
[
j
][
2
]
-=
delz
*
fpair
;
}
if
(
eflag
)
{
if
(
rsq
<
cut_coulsq
)
{
ecoul
=
prefactor
*
erfc
;
if
(
factor_coul
<
1.0
)
ecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
}
else
ecoul
=
0.0
;
if
(
rsq
<
cut_ljsq
[
itype
][
jtype
])
{
evdwl
=
r6inv
*
(
lj3
[
itype
][
jtype
]
*
r3inv
-
lj4
[
itype
][
jtype
])
-
offset
[
itype
][
jtype
];
evdwl
*=
factor_lj
;
}
else
evdwl
=
0.0
;
}
if
(
evflag
)
ev_tally
(
i
,
j
,
nlocal
,
newton_pair
,
evdwl
,
ecoul
,
fpair
,
delx
,
dely
,
delz
);
}
}
}
if
(
vflag_fdotr
)
virial_fdotr_compute
();
}
/* ----------------------------------------------------------------------
allocate all arrays
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
allocate
()
{
allocated
=
1
;
int
n
=
atom
->
ntypes
;
memory
->
create
(
setflag
,
n
+
1
,
n
+
1
,
"pair:setflag"
);
for
(
int
i
=
1
;
i
<=
n
;
i
++
)
for
(
int
j
=
i
;
j
<=
n
;
j
++
)
setflag
[
i
][
j
]
=
0
;
memory
->
create
(
cutsq
,
n
+
1
,
n
+
1
,
"pair:cutsq"
);
memory
->
create
(
cut_lj
,
n
+
1
,
n
+
1
,
"pair:cut_lj"
);
memory
->
create
(
cut_ljsq
,
n
+
1
,
n
+
1
,
"pair:cut_ljsq"
);
memory
->
create
(
epsilon
,
n
+
1
,
n
+
1
,
"pair:epsilon"
);
memory
->
create
(
sigma
,
n
+
1
,
n
+
1
,
"pair:sigma"
);
memory
->
create
(
lj1
,
n
+
1
,
n
+
1
,
"pair:lj1"
);
memory
->
create
(
lj2
,
n
+
1
,
n
+
1
,
"pair:lj2"
);
memory
->
create
(
lj3
,
n
+
1
,
n
+
1
,
"pair:lj3"
);
memory
->
create
(
lj4
,
n
+
1
,
n
+
1
,
"pair:lj4"
);
memory
->
create
(
offset
,
n
+
1
,
n
+
1
,
"pair:offset"
);
}
/* ----------------------------------------------------------------------
global settings
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
settings
(
int
narg
,
char
**
arg
)
{
if
(
narg
<
1
||
narg
>
2
)
error
->
all
(
FLERR
,
"Illegal pair_style command"
);
cut_lj_global
=
force
->
numeric
(
arg
[
0
]);
if
(
narg
==
1
)
cut_coul
=
cut_lj_global
;
else
cut_coul
=
force
->
numeric
(
arg
[
1
]);
// reset cutoffs that have been explicitly set
if
(
allocated
)
{
int
i
,
j
;
for
(
i
=
1
;
i
<=
atom
->
ntypes
;
i
++
)
for
(
j
=
i
+
1
;
j
<=
atom
->
ntypes
;
j
++
)
if
(
setflag
[
i
][
j
])
cut_lj
[
i
][
j
]
=
cut_lj_global
;
}
}
/* ----------------------------------------------------------------------
set coeffs for one or more type pairs
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
coeff
(
int
narg
,
char
**
arg
)
{
if
(
narg
<
4
||
narg
>
6
)
error
->
all
(
FLERR
,
"Incorrect args for pair coefficients"
);
if
(
!
allocated
)
allocate
();
int
ilo
,
ihi
,
jlo
,
jhi
;
force
->
bounds
(
arg
[
0
],
atom
->
ntypes
,
ilo
,
ihi
);
force
->
bounds
(
arg
[
1
],
atom
->
ntypes
,
jlo
,
jhi
);
double
epsilon_one
=
force
->
numeric
(
arg
[
2
]);
double
sigma_one
=
force
->
numeric
(
arg
[
3
]);
double
cut_lj_one
=
cut_lj_global
;
if
(
narg
==
5
)
cut_lj_one
=
force
->
numeric
(
arg
[
4
]);
int
count
=
0
;
for
(
int
i
=
ilo
;
i
<=
ihi
;
i
++
)
{
for
(
int
j
=
MAX
(
jlo
,
i
);
j
<=
jhi
;
j
++
)
{
epsilon
[
i
][
j
]
=
epsilon_one
;
sigma
[
i
][
j
]
=
sigma_one
;
cut_lj
[
i
][
j
]
=
cut_lj_one
;
setflag
[
i
][
j
]
=
1
;
count
++
;
}
}
if
(
count
==
0
)
error
->
all
(
FLERR
,
"Incorrect args for pair coefficients"
);
}
/* ----------------------------------------------------------------------
init specific to this pair style
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
init_style
()
{
if
(
!
atom
->
q_flag
)
error
->
all
(
FLERR
,
"Pair style lj/class2/coul/long requires atom attribute q"
);
neighbor
->
request
(
this
);
cut_coulsq
=
cut_coul
*
cut_coul
;
// insure use of KSpace long-range solver, set g_ewald
if
(
force
->
kspace
==
NULL
)
error
->
all
(
FLERR
,
"Pair style is incompatible with KSpace style"
);
g_ewald
=
force
->
kspace
->
g_ewald
;
}
/* ----------------------------------------------------------------------
init for one type pair i,j and corresponding j,i
------------------------------------------------------------------------- */
double
PairLJClass2CoulLong
::
init_one
(
int
i
,
int
j
)
{
// always mix epsilon,sigma via sixthpower rules
// mix distance via user-defined rule
if
(
setflag
[
i
][
j
]
==
0
)
{
epsilon
[
i
][
j
]
=
2.0
*
sqrt
(
epsilon
[
i
][
i
]
*
epsilon
[
j
][
j
])
*
pow
(
sigma
[
i
][
i
],
3.0
)
*
pow
(
sigma
[
j
][
j
],
3.0
)
/
(
pow
(
sigma
[
i
][
i
],
6.0
)
+
pow
(
sigma
[
j
][
j
],
6.0
));
sigma
[
i
][
j
]
=
pow
((
0.5
*
(
pow
(
sigma
[
i
][
i
],
6.0
)
+
pow
(
sigma
[
j
][
j
],
6.0
))),
1.0
/
6.0
);
cut_lj
[
i
][
j
]
=
mix_distance
(
cut_lj
[
i
][
i
],
cut_lj
[
j
][
j
]);
}
double
cut
=
MAX
(
cut_lj
[
i
][
j
],
cut_coul
);
cut_ljsq
[
i
][
j
]
=
cut_lj
[
i
][
j
]
*
cut_lj
[
i
][
j
];
lj1
[
i
][
j
]
=
18.0
*
epsilon
[
i
][
j
]
*
pow
(
sigma
[
i
][
j
],
9.0
);
lj2
[
i
][
j
]
=
18.0
*
epsilon
[
i
][
j
]
*
pow
(
sigma
[
i
][
j
],
6.0
);
lj3
[
i
][
j
]
=
2.0
*
epsilon
[
i
][
j
]
*
pow
(
sigma
[
i
][
j
],
9.0
);
lj4
[
i
][
j
]
=
3.0
*
epsilon
[
i
][
j
]
*
pow
(
sigma
[
i
][
j
],
6.0
);
if
(
offset_flag
)
{
double
ratio
=
sigma
[
i
][
j
]
/
cut_lj
[
i
][
j
];
offset
[
i
][
j
]
=
epsilon
[
i
][
j
]
*
(
2.0
*
pow
(
ratio
,
9.0
)
-
3.0
*
pow
(
ratio
,
6.0
));
}
else
offset
[
i
][
j
]
=
0.0
;
cut_ljsq
[
j
][
i
]
=
cut_ljsq
[
i
][
j
];
lj1
[
j
][
i
]
=
lj1
[
i
][
j
];
lj2
[
j
][
i
]
=
lj2
[
i
][
j
];
lj3
[
j
][
i
]
=
lj3
[
i
][
j
];
lj4
[
j
][
i
]
=
lj4
[
i
][
j
];
offset
[
j
][
i
]
=
offset
[
i
][
j
];
// compute I,J contribution to long-range tail correction
// count total # of atoms of type I and J via Allreduce
if
(
tail_flag
)
{
int
*
type
=
atom
->
type
;
int
nlocal
=
atom
->
nlocal
;
double
count
[
2
],
all
[
2
];
count
[
0
]
=
count
[
1
]
=
0.0
;
for
(
int
k
=
0
;
k
<
nlocal
;
k
++
)
{
if
(
type
[
k
]
==
i
)
count
[
0
]
+=
1.0
;
if
(
type
[
k
]
==
j
)
count
[
1
]
+=
1.0
;
}
MPI_Allreduce
(
count
,
all
,
2
,
MPI_DOUBLE
,
MPI_SUM
,
world
);
double
sig3
=
sigma
[
i
][
j
]
*
sigma
[
i
][
j
]
*
sigma
[
i
][
j
];
double
sig6
=
sig3
*
sig3
;
double
rc3
=
cut_lj
[
i
][
j
]
*
cut_lj
[
i
][
j
]
*
cut_lj
[
i
][
j
];
double
rc6
=
rc3
*
rc3
;
etail_ij
=
2.0
*
MY_PI
*
all
[
0
]
*
all
[
1
]
*
epsilon
[
i
][
j
]
*
sig6
*
(
sig3
-
3.0
*
rc3
)
/
(
3.0
*
rc6
);
ptail_ij
=
2.0
*
MY_PI
*
all
[
0
]
*
all
[
1
]
*
epsilon
[
i
][
j
]
*
sig6
*
(
sig3
-
2.0
*
rc3
)
/
rc6
;
}
return
cut
;
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
write_restart
(
FILE
*
fp
)
{
write_restart_settings
(
fp
);
int
i
,
j
;
for
(
i
=
1
;
i
<=
atom
->
ntypes
;
i
++
)
for
(
j
=
i
;
j
<=
atom
->
ntypes
;
j
++
)
{
fwrite
(
&
setflag
[
i
][
j
],
sizeof
(
int
),
1
,
fp
);
if
(
setflag
[
i
][
j
])
{
fwrite
(
&
epsilon
[
i
][
j
],
sizeof
(
double
),
1
,
fp
);
fwrite
(
&
sigma
[
i
][
j
],
sizeof
(
double
),
1
,
fp
);
fwrite
(
&
cut_lj
[
i
][
j
],
sizeof
(
double
),
1
,
fp
);
}
}
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
read_restart
(
FILE
*
fp
)
{
read_restart_settings
(
fp
);
allocate
();
int
i
,
j
;
int
me
=
comm
->
me
;
for
(
i
=
1
;
i
<=
atom
->
ntypes
;
i
++
)
for
(
j
=
i
;
j
<=
atom
->
ntypes
;
j
++
)
{
if
(
me
==
0
)
fread
(
&
setflag
[
i
][
j
],
sizeof
(
int
),
1
,
fp
);
MPI_Bcast
(
&
setflag
[
i
][
j
],
1
,
MPI_INT
,
0
,
world
);
if
(
setflag
[
i
][
j
])
{
if
(
me
==
0
)
{
fread
(
&
epsilon
[
i
][
j
],
sizeof
(
double
),
1
,
fp
);
fread
(
&
sigma
[
i
][
j
],
sizeof
(
double
),
1
,
fp
);
fread
(
&
cut_lj
[
i
][
j
],
sizeof
(
double
),
1
,
fp
);
}
MPI_Bcast
(
&
epsilon
[
i
][
j
],
1
,
MPI_DOUBLE
,
0
,
world
);
MPI_Bcast
(
&
sigma
[
i
][
j
],
1
,
MPI_DOUBLE
,
0
,
world
);
MPI_Bcast
(
&
cut_lj
[
i
][
j
],
1
,
MPI_DOUBLE
,
0
,
world
);
}
}
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
write_restart_settings
(
FILE
*
fp
)
{
fwrite
(
&
cut_lj_global
,
sizeof
(
double
),
1
,
fp
);
fwrite
(
&
cut_coul
,
sizeof
(
double
),
1
,
fp
);
fwrite
(
&
offset_flag
,
sizeof
(
int
),
1
,
fp
);
fwrite
(
&
mix_flag
,
sizeof
(
int
),
1
,
fp
);
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void
PairLJClass2CoulLong
::
read_restart_settings
(
FILE
*
fp
)
{
if
(
comm
->
me
==
0
)
{
fread
(
&
cut_lj_global
,
sizeof
(
double
),
1
,
fp
);
fread
(
&
cut_coul
,
sizeof
(
double
),
1
,
fp
);
fread
(
&
offset_flag
,
sizeof
(
int
),
1
,
fp
);
fread
(
&
mix_flag
,
sizeof
(
int
),
1
,
fp
);
}
MPI_Bcast
(
&
cut_lj_global
,
1
,
MPI_DOUBLE
,
0
,
world
);
MPI_Bcast
(
&
cut_coul
,
1
,
MPI_DOUBLE
,
0
,
world
);
MPI_Bcast
(
&
offset_flag
,
1
,
MPI_INT
,
0
,
world
);
MPI_Bcast
(
&
mix_flag
,
1
,
MPI_INT
,
0
,
world
);
}
/* ---------------------------------------------------------------------- */
double
PairLJClass2CoulLong
::
single
(
int
i
,
int
j
,
int
itype
,
int
jtype
,
double
rsq
,
double
factor_coul
,
double
factor_lj
,
double
&
fforce
)
{
double
r2inv
,
r
,
rinv
,
r3inv
,
r6inv
,
grij
,
expm2
,
t
,
erfc
,
prefactor
;
double
forcecoul
,
forcelj
,
phicoul
,
philj
;
r2inv
=
1.0
/
rsq
;
if
(
rsq
<
cut_coulsq
)
{
r
=
sqrt
(
rsq
);
grij
=
g_ewald
*
r
;
expm2
=
exp
(
-
grij
*
grij
);
t
=
1.0
/
(
1.0
+
EWALD_P
*
grij
);
erfc
=
t
*
(
A1
+
t
*
(
A2
+
t
*
(
A3
+
t
*
(
A4
+
t
*
A5
))))
*
expm2
;
prefactor
=
force
->
qqrd2e
*
atom
->
q
[
i
]
*
atom
->
q
[
j
]
/
r
;
forcecoul
=
prefactor
*
(
erfc
+
EWALD_F
*
grij
*
expm2
);
if
(
factor_coul
<
1.0
)
forcecoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
}
else
forcecoul
=
0.0
;
if
(
rsq
<
cut_ljsq
[
itype
][
jtype
])
{
rinv
=
sqrt
(
r2inv
);
r3inv
=
r2inv
*
rinv
;
r6inv
=
r3inv
*
r3inv
;
forcelj
=
r6inv
*
(
lj1
[
itype
][
jtype
]
*
r3inv
-
lj2
[
itype
][
jtype
]);
}
else
forcelj
=
0.0
;
fforce
=
(
forcecoul
+
factor_lj
*
forcelj
)
*
r2inv
;
double
eng
=
0.0
;
if
(
rsq
<
cut_coulsq
)
{
phicoul
=
prefactor
*
erfc
;
if
(
factor_coul
<
1.0
)
phicoul
-=
(
1.0
-
factor_coul
)
*
prefactor
;
eng
+=
phicoul
;
}
if
(
rsq
<
cut_ljsq
[
itype
][
jtype
])
{
philj
=
r6inv
*
(
lj3
[
itype
][
jtype
]
*
r3inv
-
lj4
[
itype
][
jtype
])
-
offset
[
itype
][
jtype
];
eng
+=
factor_lj
*
philj
;
}
return
eng
;
}
/* ---------------------------------------------------------------------- */
void
*
PairLJClass2CoulLong
::
extract
(
const
char
*
str
,
int
&
dim
)
{
dim
=
0
;
if
(
strcmp
(
str
,
"cut_coul"
)
==
0
)
return
(
void
*
)
&
cut_coul
;
return
NULL
;
}
Event Timeline
Log In to Comment