-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathbond_fene_expand.cpp
275 lines (222 loc) · 8.03 KB
/
bond_fene_expand.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
/* ----------------------------------------------------------------------
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 "stdlib.h"
#include "bond_fene_expand.h"
#include "atom.h"
#include "neighbor.h"
#include "domain.h"
#include "comm.h"
#include "update.h"
#include "force.h"
#include "memory.h"
#include "error.h"
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
BondFENEExpand::BondFENEExpand(LAMMPS *lmp) : Bond(lmp)
{
TWO_1_3 = pow(2.0,(1.0/3.0));
}
/* ---------------------------------------------------------------------- */
BondFENEExpand::~BondFENEExpand()
{
if (allocated) {
memory->sfree(setflag);
memory->sfree(k);
memory->sfree(r0);
memory->sfree(epsilon);
memory->sfree(sigma);
memory->sfree(shift);
}
}
/* ---------------------------------------------------------------------- */
void BondFENEExpand::compute(int eflag, int vflag)
{
int i1,i2,n,type;
double delx,dely,delz,ebond,fbond;
double rsq,r0sq,rlogarg,sr2,sr6;
double r,rshift,rshiftsq;
ebond = 0.0;
if (eflag || vflag) ev_setup(eflag,vflag);
else evflag = 0;
double **x = atom->x;
double **f = atom->f;
int **bondlist = neighbor->bondlist;
int nbondlist = neighbor->nbondlist;
int nlocal = atom->nlocal;
int newton_bond = force->newton_bond;
for (n = 0; n < nbondlist; n++) {
i1 = bondlist[n][0];
i2 = bondlist[n][1];
type = bondlist[n][2];
delx = x[i1][0] - x[i2][0];
dely = x[i1][1] - x[i2][1];
delz = x[i1][2] - x[i2][2];
domain->minimum_image(delx,dely,delz);
// force from log term
rsq = delx*delx + dely*dely + delz*delz;
r = sqrt(rsq);
rshift = r - shift[type];
rshiftsq = rshift*rshift;
r0sq = r0[type] * r0[type];
rlogarg = 1.0 - rshiftsq/r0sq;
// if r -> r0, then rlogarg < 0.0 which is an error
// issue a warning and reset rlogarg = epsilon
// if r > 2*r0 something serious is wrong, abort
if (rlogarg < 0.1) {
char str[128];
sprintf(str,"FENE bond too long: %d %d %d %g",
update->ntimestep,atom->tag[i1],atom->tag[i2],sqrt(rsq));
error->warning(str);
if (rlogarg <= -3.0) error->one("Bad FENE bond");
rlogarg = 0.1;
}
fbond = -k[type]*rshift/rlogarg/r;
// force from LJ term
if (rshiftsq < TWO_1_3*sigma[type]*sigma[type]) {
sr2 = sigma[type]*sigma[type]/rshiftsq;
sr6 = sr2*sr2*sr2;
fbond += 48.0*epsilon[type]*sr6*(sr6-0.5)/rshift/r;
}
// energy
if (eflag) {
ebond = -0.5 * k[type]*r0sq*log(rlogarg);
if (rshiftsq < TWO_1_3*sigma[type]*sigma[type])
ebond += 4.0*epsilon[type]*sr6*(sr6-1.0) + epsilon[type];
}
// apply force to each of 2 atoms
if (newton_bond || i1 < nlocal) {
f[i1][0] += delx*fbond;
f[i1][1] += dely*fbond;
f[i1][2] += delz*fbond;
}
if (newton_bond || i2 < nlocal) {
f[i2][0] -= delx*fbond;
f[i2][1] -= dely*fbond;
f[i2][2] -= delz*fbond;
}
if (evflag) ev_tally(i1,i2,nlocal,newton_bond,ebond,fbond,delx,dely,delz);
}
}
/* ---------------------------------------------------------------------- */
void BondFENEExpand::allocate()
{
allocated = 1;
int n = atom->nbondtypes;
k = (double *) memory->smalloc((n+1)*sizeof(double),"bond:k");
r0 = (double *) memory->smalloc((n+1)*sizeof(double),"bond:r0");
epsilon = (double *) memory->smalloc((n+1)*sizeof(double),"bond:epsilon");
sigma = (double *) memory->smalloc((n+1)*sizeof(double),"bond:sigma");
shift = (double *) memory->smalloc((n+1)*sizeof(double),"bond:shift");
setflag = (int *) memory->smalloc((n+1)*sizeof(int),"bond:setflag");
for (int i = 1; i <= n; i++) setflag[i] = 0;
}
/* ----------------------------------------------------------------------
set coeffs for one type
------------------------------------------------------------------------- */
void BondFENEExpand::coeff(int narg, char **arg)
{
if (narg != 6) error->all("Incorrect args for bond coefficients");
if (!allocated) allocate();
int ilo,ihi;
force->bounds(arg[0],atom->nbondtypes,ilo,ihi);
double k_one = force->numeric(arg[1]);
double r0_one = force->numeric(arg[2]);
double epsilon_one = force->numeric(arg[3]);
double sigma_one = force->numeric(arg[4]);
double shift_one = force->numeric(arg[5]);
int count = 0;
for (int i = ilo; i <= ihi; i++) {
k[i] = k_one;
r0[i] = r0_one;
epsilon[i] = epsilon_one;
sigma[i] = sigma_one;
shift[i] = shift_one;
setflag[i] = 1;
count++;
}
if (count == 0) error->all("Incorrect args for bond coefficients");
}
/* ----------------------------------------------------------------------
check if special_bond settings are valid
------------------------------------------------------------------------- */
void BondFENEExpand::init_style()
{
// special bonds should be 0 1 1
if (force->special_lj[1] != 0.0 || force->special_lj[2] != 1.0 ||
force->special_lj[3] != 1.0) {
if (comm->me == 0)
error->warning("Use special bonds = 0,1,1 with bond style fene/expand");
}
}
/* ---------------------------------------------------------------------- */
double BondFENEExpand::equilibrium_distance(int i)
{
return 0.97*sigma[i] + shift[i];
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void BondFENEExpand::write_restart(FILE *fp)
{
fwrite(&k[1],sizeof(double),atom->nbondtypes,fp);
fwrite(&r0[1],sizeof(double),atom->nbondtypes,fp);
fwrite(&epsilon[1],sizeof(double),atom->nbondtypes,fp);
fwrite(&sigma[1],sizeof(double),atom->nbondtypes,fp);
fwrite(&shift[1],sizeof(double),atom->nbondtypes,fp);
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void BondFENEExpand::read_restart(FILE *fp)
{
allocate();
if (comm->me == 0) {
fread(&k[1],sizeof(double),atom->nbondtypes,fp);
fread(&r0[1],sizeof(double),atom->nbondtypes,fp);
fread(&epsilon[1],sizeof(double),atom->nbondtypes,fp);
fread(&sigma[1],sizeof(double),atom->nbondtypes,fp);
fread(&shift[1],sizeof(double),atom->nbondtypes,fp);
}
MPI_Bcast(&k[1],atom->nbondtypes,MPI_DOUBLE,0,world);
MPI_Bcast(&r0[1],atom->nbondtypes,MPI_DOUBLE,0,world);
MPI_Bcast(&epsilon[1],atom->nbondtypes,MPI_DOUBLE,0,world);
MPI_Bcast(&sigma[1],atom->nbondtypes,MPI_DOUBLE,0,world);
MPI_Bcast(&shift[1],atom->nbondtypes,MPI_DOUBLE,0,world);
for (int i = 1; i <= atom->nbondtypes; i++) setflag[i] = 1;
}
/* ---------------------------------------------------------------------- */
double BondFENEExpand::single(int type, double rsq, int i, int j)
{
double r = sqrt(rsq);
double rshift = r - shift[type];
double rshiftsq = rshift*rshift;
double r0sq = r0[type] * r0[type];
double rlogarg = 1.0 - rshiftsq/r0sq;
// if r -> r0, then rlogarg < 0.0 which is an error
// issue a warning and reset rlogarg = epsilon
// if r > 2*r0 something serious is wrong, abort
if (rlogarg < 0.1) {
char str[128];
sprintf(str,"FENE bond too long: %d %g",update->ntimestep,sqrt(rsq));
error->warning(str);
if (rlogarg <= -3.0) error->one("Bad FENE bond");
rlogarg = 0.1;
}
double eng = -0.5 * k[type]*r0sq*log(rlogarg);
if (rshiftsq < TWO_1_3*sigma[type]*sigma[type]) {
double sr2,sr6;
sr2 = sigma[type]*sigma[type]/rshiftsq;
sr6 = sr2*sr2*sr2;
eng += 4.0*epsilon[type]*sr6*(sr6-1.0) + epsilon[type];
}
return eng;
}