forked from lammps/lammps
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpair_tersoff_mod.cpp
338 lines (269 loc) · 11.5 KB
/
pair_tersoff_mod.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
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, [email protected]
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.
------------------------------------------------------------------------- */
/* ----------------------------------------------------------------------
Contributing author: Aidan Thompson (SNL) - original Tersoff implementation
Vitaly Dozhdikov (JIHT of RAS) - MOD addition
------------------------------------------------------------------------- */
#include "pair_tersoff_mod.h"
#include <mpi.h>
#include <cmath>
#include <cstdlib>
#include <cstring>
#include "atom.h"
#include "force.h"
#include "comm.h"
#include "math_const.h"
#include "math_special.h"
#include "memory.h"
#include "error.h"
#include "utils.h"
#include "tokenizer.h"
#include "potential_file_reader.h"
using namespace LAMMPS_NS;
using namespace MathConst;
using namespace MathSpecial;
#define DELTA 4
/* ---------------------------------------------------------------------- */
PairTersoffMOD::PairTersoffMOD(LAMMPS *lmp) : PairTersoff(lmp) {}
/* ---------------------------------------------------------------------- */
void PairTersoffMOD::read_file(char *file)
{
memory->sfree(params);
params = nullptr;
nparams = maxparam = 0;
// open file on proc 0
if (comm->me == 0) {
PotentialFileReader reader(lmp, file, "tersoff/mod", unit_convert_flag);
char * line;
// transparently convert units for supported conversions
int unit_convert = reader.get_unit_convert();
double conversion_factor = utils::get_conversion_factor(utils::ENERGY,
unit_convert);
while((line = reader.next_line(NPARAMS_PER_LINE))) {
try {
ValueTokenizer values(line);
std::string iname = values.next_string();
std::string jname = values.next_string();
std::string kname = values.next_string();
// ielement,jelement,kelement = 1st args
// if all 3 args are in element list, then parse this line
// else skip to next entry in file
int ielement, jelement, kelement;
for (ielement = 0; ielement < nelements; ielement++)
if (iname == elements[ielement]) break;
if (ielement == nelements) continue;
for (jelement = 0; jelement < nelements; jelement++)
if (jname == elements[jelement]) break;
if (jelement == nelements) continue;
for (kelement = 0; kelement < nelements; kelement++)
if (kname == elements[kelement]) break;
if (kelement == nelements) continue;
// load up parameter settings and error check their values
if (nparams == maxparam) {
maxparam += DELTA;
params = (Param *) memory->srealloc(params,maxparam*sizeof(Param),
"pair:params");
// make certain all addional allocated storage is initialized
// to avoid false positives when checking with valgrind
memset(params + nparams, 0, DELTA*sizeof(Param));
}
params[nparams].ielement = ielement;
params[nparams].jelement = jelement;
params[nparams].kelement = kelement;
params[nparams].powerm = values.next_double();
params[nparams].lam3 = values.next_double();
params[nparams].h = values.next_double();
params[nparams].powern = values.next_double();
params[nparams].beta = values.next_double();
params[nparams].lam2 = values.next_double();
params[nparams].bigb = values.next_double();
params[nparams].bigr = values.next_double();
params[nparams].bigd = values.next_double();
params[nparams].lam1 = values.next_double();
params[nparams].biga = values.next_double();
params[nparams].powern_del = values.next_double();
params[nparams].c1 = values.next_double();
params[nparams].c2 = values.next_double();
params[nparams].c3 = values.next_double();
params[nparams].c4 = values.next_double();
params[nparams].c5 = values.next_double();
params[nparams].powermint = int(params[nparams].powerm);
if (unit_convert) {
params[nparams].biga *= conversion_factor;
params[nparams].bigb *= conversion_factor;
}
} catch (TokenizerException & e) {
error->one(FLERR, e.what());
}
// currently only allow m exponent of 1 or 3
if (params[nparams].powern < 0.0 ||
params[nparams].beta < 0.0 ||
params[nparams].lam2 < 0.0 ||
params[nparams].bigb < 0.0 ||
params[nparams].bigr < 0.0 ||
params[nparams].bigd < 0.0 ||
params[nparams].bigd > params[nparams].bigr ||
params[nparams].lam1 < 0.0 ||
params[nparams].biga < 0.0 ||
params[nparams].powerm - params[nparams].powermint != 0.0 ||
(params[nparams].powermint != 3 &&
params[nparams].powermint != 1)
)
error->one(FLERR,"Illegal Tersoff parameter");
nparams++;
}
}
MPI_Bcast(&nparams, 1, MPI_INT, 0, world);
MPI_Bcast(&maxparam, 1, MPI_INT, 0, world);
if(comm->me != 0) {
params = (Param *) memory->srealloc(params,maxparam*sizeof(Param), "pair:params");
}
MPI_Bcast(params, maxparam*sizeof(Param), MPI_BYTE, 0, world);
}
/* ---------------------------------------------------------------------- */
void PairTersoffMOD::setup_params()
{
int i,j,k,m,n;
// set elem2param for all element triplet combinations
// must be a single exact match to lines read from file
// do not allow for ACB in place of ABC
memory->destroy(elem2param);
memory->create(elem2param,nelements,nelements,nelements,"pair:elem2param");
for (i = 0; i < nelements; i++)
for (j = 0; j < nelements; j++)
for (k = 0; k < nelements; k++) {
n = -1;
for (m = 0; m < nparams; m++) {
if (i == params[m].ielement && j == params[m].jelement &&
k == params[m].kelement) {
if (n >= 0) error->all(FLERR,"Potential file has duplicate entry");
n = m;
}
}
if (n < 0) error->all(FLERR,"Potential file is missing an entry");
elem2param[i][j][k] = n;
}
// compute parameter values derived from inputs
for (m = 0; m < nparams; m++) {
params[m].cut = params[m].bigr + params[m].bigd;
params[m].cutsq = params[m].cut*params[m].cut;
if (params[m].powern > 0.0) {
params[m].ca1 = pow(2.0*params[m].powern_del*1.0e-16,-1.0/params[m].powern);
params[m].ca4 = 1.0/params[m].ca1;
} else params[m].ca1 = params[m].ca4 = 0.0;
}
// set cutmax to max of all params
cutmax = 0.0;
for (m = 0; m < nparams; m++)
if (params[m].cut > cutmax) cutmax = params[m].cut;
}
/* ---------------------------------------------------------------------- */
double PairTersoffMOD::zeta(Param *param, double rsqij, double rsqik,
double *delrij, double *delrik)
{
double rij,rik,costheta,arg,ex_delr;
rij = sqrt(rsqij);
rik = sqrt(rsqik);
costheta = (delrij[0]*delrik[0] + delrij[1]*delrik[1] +
delrij[2]*delrik[2]) / (rij*rik);
if (param->powermint == 3) arg = cube(param->lam3 * (rij-rik));
else arg = param->lam3 * (rij-rik);
if (arg > 69.0776) ex_delr = 1.e30;
else if (arg < -69.0776) ex_delr = 0.0;
else ex_delr = exp(arg);
return ters_fc(rik,param) * ters_gijk_mod(costheta,param) * ex_delr;
}
/* ---------------------------------------------------------------------- */
double PairTersoffMOD::ters_fc(double r, Param *param)
{
double ters_R = param->bigr;
double ters_D = param->bigd;
if (r < ters_R-ters_D) return 1.0;
if (r > ters_R+ters_D) return 0.0;
return 0.5*(1.0 - 1.125*sin(MY_PI2*(r - ters_R)/ters_D) -
0.125*sin(3*MY_PI2*(r - ters_R)/ters_D));
}
/* ---------------------------------------------------------------------- */
double PairTersoffMOD::ters_fc_d(double r, Param *param)
{
double ters_R = param->bigr;
double ters_D = param->bigd;
if (r < ters_R-ters_D) return 0.0;
if (r > ters_R+ters_D) return 0.0;
return -(0.375*MY_PI4/ters_D) * (3*cos(MY_PI2*(r - ters_R)/ters_D) +
cos(3*MY_PI2*(r - ters_R)/ters_D));
}
/* ---------------------------------------------------------------------- */
double PairTersoffMOD::ters_bij(double zeta, Param *param)
{
double tmp = param->beta * zeta;
if (tmp > param->ca1) return pow(tmp, -param->powern/(2.0*param->powern_del));
if (tmp < param->ca4) return 1.0;
return pow(1.0 + pow(tmp,param->powern), -1.0/(2.0*param->powern_del));
}
/* ---------------------------------------------------------------------- */
double PairTersoffMOD::ters_bij_d(double zeta, Param *param)
{
double tmp = param->beta * zeta;
if (tmp > param->ca1) return -0.5*(param->powern/param->powern_del)*
pow(tmp,-0.5*(param->powern/param->powern_del)) / zeta;
if (tmp < param->ca4) return 0.0;
double tmp_n = pow(tmp,param->powern);
return -0.5 *(param->powern/param->powern_del)*
pow(1.0+tmp_n, -1.0-(1.0/(2.0*param->powern_del)))*tmp_n / zeta;
}
/* ---------------------------------------------------------------------- */
void PairTersoffMOD::ters_zetaterm_d(double prefactor,
double *rij_hat, double rij,
double *rik_hat, double rik,
double *dri, double *drj, double *drk,
Param *param)
{
double gijk,gijk_d,ex_delr,ex_delr_d,fc,dfc,cos_theta,tmp;
double dcosdri[3],dcosdrj[3],dcosdrk[3];
fc = ters_fc(rik,param);
dfc = ters_fc_d(rik,param);
if (param->powermint == 3) tmp = cube(param->lam3 * (rij-rik));
else tmp = param->lam3 * (rij-rik);
if (tmp > 69.0776) ex_delr = 1.e30;
else if (tmp < -69.0776) ex_delr = 0.0;
else ex_delr = exp(tmp);
if (param->powermint == 3)
ex_delr_d = 3.0*cube(param->lam3) * square(rij-rik)*ex_delr;
else ex_delr_d = param->lam3 * ex_delr;
cos_theta = vec3_dot(rij_hat,rik_hat);
gijk = ters_gijk_mod(cos_theta,param);
gijk_d = ters_gijk_d_mod(cos_theta,param);
costheta_d(rij_hat,rij,rik_hat,rik,dcosdri,dcosdrj,dcosdrk);
// compute the derivative wrt Ri
// dri = -dfc*gijk*ex_delr*rik_hat;
// dri += fc*gijk_d*ex_delr*dcosdri;
// dri += fc*gijk*ex_delr_d*(rik_hat - rij_hat);
vec3_scale(-dfc*gijk*ex_delr,rik_hat,dri);
vec3_scaleadd(fc*gijk_d*ex_delr,dcosdri,dri,dri);
vec3_scaleadd(fc*gijk*ex_delr_d,rik_hat,dri,dri);
vec3_scaleadd(-fc*gijk*ex_delr_d,rij_hat,dri,dri);
vec3_scale(prefactor,dri,dri);
// compute the derivative wrt Rj
// drj = fc*gijk_d*ex_delr*dcosdrj;
// drj += fc*gijk*ex_delr_d*rij_hat;
vec3_scale(fc*gijk_d*ex_delr,dcosdrj,drj);
vec3_scaleadd(fc*gijk*ex_delr_d,rij_hat,drj,drj);
vec3_scale(prefactor,drj,drj);
// compute the derivative wrt Rk
// drk = dfc*gijk*ex_delr*rik_hat;
// drk += fc*gijk_d*ex_delr*dcosdrk;
// drk += -fc*gijk*ex_delr_d*rik_hat;
vec3_scale(dfc*gijk*ex_delr,rik_hat,drk);
vec3_scaleadd(fc*gijk_d*ex_delr,dcosdrk,drk,drk);
vec3_scaleadd(-fc*gijk*ex_delr_d,rik_hat,drk,drk);
vec3_scale(prefactor,drk,drk);
}