forked from lammps/lammps
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathreaxff_lookup.cpp
291 lines (246 loc) · 10.6 KB
/
reaxff_lookup.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
// clang-format off
/*----------------------------------------------------------------------
PuReMD - Purdue ReaxFF Molecular Dynamics Program
Copyright (2010) Purdue University
Hasan Metin Aktulga, [email protected]
Joseph Fogarty, [email protected]
Sagar Pandit, [email protected]
Ananth Y Grama, [email protected]
Please cite the related publication:
H. M. Aktulga, J. C. Fogarty, S. A. Pandit, A. Y. Grama,
"Parallel Reactive Molecular Dynamics: Numerical Methods and
Algorithmic Techniques", Parallel Computing, in press.
This program is free software; you can redistribute it and/or
modify it under the terms of the GNU General Public License as
published by the Free Software Foundation; either version 2 of
the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
See the GNU General Public License for more details:
<https://www.gnu.org/licenses/>.
----------------------------------------------------------------------*/
#include "reaxff_api.h"
#include <mpi.h>
#include <cstdlib>
namespace ReaxFF {
static void Tridiagonal_Solve(const double *a, const double *b,
double *c, double *d, double *x, unsigned int n) {
int i;
double id;
c[0] /= b[0]; /* Division by zero risk. */
d[0] /= b[0]; /* Division by zero would imply a singular matrix. */
for (i = 1; i < (int)n; i++) {
id = (b[i] - c[i-1] * a[i]); /* Division by zero risk. */
c[i] /= id; /* Last value calculated is redundant. */
d[i] = (d[i] - d[i-1] * a[i])/id;
}
x[n - 1] = d[n - 1];
for (i = n - 2; i >= 0; i--)
x[i] = d[i] - c[i] * x[i + 1];
}
void Natural_Cubic_Spline(LAMMPS_NS::Error* error_ptr, const double *h, const double *f,
cubic_spline_coef *coef, unsigned int n)
{
int i;
double *a, *b, *c, *d, *v;
/* allocate space for the linear system */
a = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
b = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
c = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
d = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
v = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
/* build the linear system */
a[0] = a[1] = a[n-1] = 0;
for (i = 2; i < (int)n-1; ++i)
a[i] = h[i-1];
b[0] = b[n-1] = 0;
for (i = 1; i < (int)n-1; ++i)
b[i] = 2 * (h[i-1] + h[i]);
c[0] = c[n-2] = c[n-1] = 0;
for (i = 1; i < (int)n-2; ++i)
c[i] = h[i];
d[0] = d[n-1] = 0;
for (i = 1; i < (int)n-1; ++i)
d[i] = 6 * ((f[i+1]-f[i])/h[i] - (f[i]-f[i-1])/h[i-1]);
v[0] = 0;
v[n-1] = 0;
Tridiagonal_Solve(&(a[1]), &(b[1]), &(c[1]), &(d[1]), &(v[1]), n-2);
for (i = 1; i < (int)n; ++i) {
coef[i-1].d = (v[i] - v[i-1]) / (6*h[i-1]);
coef[i-1].c = v[i]/2;
coef[i-1].b = (f[i]-f[i-1])/h[i-1] + h[i-1]*(2*v[i] + v[i-1])/6;
coef[i-1].a = f[i];
}
sfree(error_ptr, a, "cubic_spline:a");
sfree(error_ptr, b, "cubic_spline:b");
sfree(error_ptr, c, "cubic_spline:c");
sfree(error_ptr, d, "cubic_spline:d");
sfree(error_ptr, v, "cubic_spline:v");
}
void Complete_Cubic_Spline(LAMMPS_NS::Error* error_ptr, const double *h,
const double *f, double v0, double vlast,
cubic_spline_coef *coef, unsigned int n)
{
int i;
double *a, *b, *c, *d, *v;
/* allocate space for the linear system */
a = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
b = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
c = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
d = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
v = (double*) smalloc(error_ptr, n * sizeof(double), "cubic_spline:a");
/* build the linear system */
a[0] = 0;
for (i = 1; i < (int)n; ++i)
a[i] = h[i-1];
b[0] = 2*h[0];
for (i = 1; i < (int)n; ++i)
b[i] = 2 * (h[i-1] + h[i]);
c[n-1] = 0;
for (i = 0; i < (int)n-1; ++i)
c[i] = h[i];
d[0] = 6 * (f[1]-f[0])/h[0] - 6 * v0;
d[n-1] = 6 * vlast - 6 * (f[n-1]-f[n-2]/h[n-2]);
for (i = 1; i < (int)n-1; ++i)
d[i] = 6 * ((f[i+1]-f[i])/h[i] - (f[i]-f[i-1])/h[i-1]);
Tridiagonal_Solve(&(a[0]), &(b[0]), &(c[0]), &(d[0]), &(v[0]), n);
for (i = 1; i < (int)n; ++i) {
coef[i-1].d = (v[i] - v[i-1]) / (6*h[i-1]);
coef[i-1].c = v[i]/2;
coef[i-1].b = (f[i]-f[i-1])/h[i-1] + h[i-1]*(2*v[i] + v[i-1])/6;
coef[i-1].a = f[i];
}
sfree(error_ptr, a, "cubic_spline:a");
sfree(error_ptr, b, "cubic_spline:b");
sfree(error_ptr, c, "cubic_spline:c");
sfree(error_ptr, d, "cubic_spline:d");
sfree(error_ptr, v, "cubic_spline:v");
}
void Init_Lookup_Tables(reax_system *system, control_params *control,
storage *workspace, MPI_Comm world)
{
int i, j, r;
double dr;
double *h, *fh, *fvdw, *fele, *fCEvd, *fCEclmb;
double v0_vdw, v0_ele, vlast_vdw, vlast_ele;
LR_lookup_table ** & LR = system->LR;
/* initializations */
v0_vdw = 0;
v0_ele = 0;
vlast_vdw = 0;
vlast_ele = 0;
const int num_atom_types = system->reax_param.num_atom_types;
int *existing_types = new int[num_atom_types];
int *aggregated = new int[num_atom_types];
dr = control->nonb_cut / control->tabulate;
h = (double*)
smalloc(system->error_ptr, (control->tabulate+2) * sizeof(double), "lookup:h");
fh = (double*)
smalloc(system->error_ptr, (control->tabulate+2) * sizeof(double), "lookup:fh");
fvdw = (double*)
smalloc(system->error_ptr, (control->tabulate+2) * sizeof(double), "lookup:fvdw");
fCEvd = (double*)
smalloc(system->error_ptr, (control->tabulate+2) * sizeof(double), "lookup:fCEvd");
fele = (double*)
smalloc(system->error_ptr, (control->tabulate+2) * sizeof(double), "lookup:fele");
fCEclmb = (double*)
smalloc(system->error_ptr, (control->tabulate+2) * sizeof(double), "lookup:fCEclmb");
LR = (LR_lookup_table**)
scalloc(system->error_ptr, num_atom_types, sizeof(LR_lookup_table*), "lookup:LR");
for (i = 0; i < num_atom_types; ++i)
LR[i] = (LR_lookup_table*)
scalloc(system->error_ptr, num_atom_types, sizeof(LR_lookup_table), "lookup:LR[i]");
for (i = 0; i < num_atom_types; ++i)
existing_types[i] = 0;
for (i = 0; i < system->n; ++i)
existing_types[system->my_atoms[i].type] = 1;
MPI_Allreduce(existing_types, aggregated, num_atom_types, MPI_INT, MPI_SUM, world);
for (i = 0; i < num_atom_types; ++i) {
if (aggregated[i]) {
for (j = i; j < num_atom_types; ++j) {
if (aggregated[j]) {
LR[i][j].xmin = 0;
LR[i][j].xmax = control->nonb_cut;
LR[i][j].n = control->tabulate + 2;
LR[i][j].dx = dr;
LR[i][j].inv_dx = control->tabulate / control->nonb_cut;
LR[i][j].y = (LR_data*)
smalloc(system->error_ptr, LR[i][j].n * sizeof(LR_data), "lookup:LR[i,j].y");
LR[i][j].H = (cubic_spline_coef*)
smalloc(system->error_ptr, LR[i][j].n*sizeof(cubic_spline_coef),"lookup:LR[i,j].H");
LR[i][j].vdW = (cubic_spline_coef*)
smalloc(system->error_ptr, LR[i][j].n*sizeof(cubic_spline_coef),"lookup:LR[i,j].vdW");
LR[i][j].CEvd = (cubic_spline_coef*)
smalloc(system->error_ptr, LR[i][j].n*sizeof(cubic_spline_coef),"lookup:LR[i,j].CEvd");
LR[i][j].ele = (cubic_spline_coef*)
smalloc(system->error_ptr, LR[i][j].n*sizeof(cubic_spline_coef),"lookup:LR[i,j].ele");
LR[i][j].CEclmb = (cubic_spline_coef*)
smalloc(system->error_ptr, LR[i][j].n*sizeof(cubic_spline_coef),
"lookup:LR[i,j].CEclmb");
for (r = 1; r <= control->tabulate; ++r) {
LR_vdW_Coulomb(system, workspace, control, i, j, r * dr, &(LR[i][j].y[r]));
h[r] = LR[i][j].dx;
fh[r] = LR[i][j].y[r].H;
fvdw[r] = LR[i][j].y[r].e_vdW;
fCEvd[r] = LR[i][j].y[r].CEvd;
fele[r] = LR[i][j].y[r].e_ele;
fCEclmb[r] = LR[i][j].y[r].CEclmb;
}
// init the start-end points
h[r] = LR[i][j].dx;
v0_vdw = LR[i][j].y[1].CEvd;
v0_ele = LR[i][j].y[1].CEclmb;
fh[r] = fh[r-1];
fvdw[r] = fvdw[r-1];
fCEvd[r] = fCEvd[r-1];
fele[r] = fele[r-1];
fCEclmb[r] = fCEclmb[r-1];
vlast_vdw = fCEvd[r-1];
vlast_ele = fele[r-1];
Natural_Cubic_Spline(control->error_ptr, &h[1], &fh[1],
&(LR[i][j].H[1]), control->tabulate+1);
Complete_Cubic_Spline(control->error_ptr, &h[1], &fvdw[1], v0_vdw, vlast_vdw,
&(LR[i][j].vdW[1]), control->tabulate+1);
Natural_Cubic_Spline(control->error_ptr, &h[1], &fCEvd[1],
&(LR[i][j].CEvd[1]), control->tabulate+1);
Complete_Cubic_Spline(control->error_ptr, &h[1], &fele[1], v0_ele, vlast_ele,
&(LR[i][j].ele[1]), control->tabulate+1);
Natural_Cubic_Spline(control->error_ptr, &h[1], &fCEclmb[1],
&(LR[i][j].CEclmb[1]), control->tabulate+1);
} else {
LR[i][j].n = 0;
}
}
}
}
free(h);
free(fh);
free(fvdw);
free(fCEvd);
free(fele);
free(fCEclmb);
delete[] existing_types;
delete[] aggregated;
}
void Deallocate_Lookup_Tables(reax_system *system)
{
int i, j;
int ntypes;
LR_lookup_table ** & LR = system->LR;
ntypes = system->reax_param.num_atom_types;
for (i = 0; i < ntypes; ++i) {
for (j = i; j < ntypes; ++j)
if (LR[i][j].n) {
sfree(system->error_ptr, LR[i][j].y, "LR[i,j].y");
sfree(system->error_ptr, LR[i][j].H, "LR[i,j].H");
sfree(system->error_ptr, LR[i][j].vdW, "LR[i,j].vdW");
sfree(system->error_ptr, LR[i][j].CEvd, "LR[i,j].CEvd");
sfree(system->error_ptr, LR[i][j].ele, "LR[i,j].ele");
sfree(system->error_ptr, LR[i][j].CEclmb, "LR[i,j].CEclmb");
}
sfree(system->error_ptr, LR[i], "LR[i]");
}
sfree(system->error_ptr, LR, "LR");
}
}