1 // clang-format off
2 /* ----------------------------------------------------------------------
3 LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
4 https://www.lammps.org/, Sandia National Laboratories
5 Steve Plimpton, sjplimp@sandia.gov
6
7 Copyright (2003) Sandia Corporation. Under the terms of Contract
8 DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
9 certain rights in this software. This software is distributed under
10 the GNU General Public License.
11
12 See the README file in the top-level LAMMPS directory.
13 ------------------------------------------------------------------------- */
14
15 /* ----------------------------------------------------------------------
16 Contributing author: Pieter J. in 't Veld (SNL)
17 Tabulation for long-range dispersion added by Wayne Mitchell (Loyola
18 University New Orleans)
19 ------------------------------------------------------------------------- */
20
21 #include "pair_lj_long_coul_long.h"
22
23 #include "atom.h"
24 #include "comm.h"
25 #include "error.h"
26 #include "force.h"
27 #include "kspace.h"
28 #include "math_extra.h"
29 #include "memory.h"
30 #include "neigh_list.h"
31 #include "neigh_request.h"
32 #include "neighbor.h"
33 #include "respa.h"
34 #include "update.h"
35
36 #include <cmath>
37 #include <cstring>
38
39 using namespace LAMMPS_NS;
40 using namespace MathExtra;
41
42 #define EWALD_F 1.12837917
43 #define EWALD_P 0.3275911
44 #define A1 0.254829592
45 #define A2 -0.284496736
46 #define A3 1.421413741
47 #define A4 -1.453152027
48 #define A5 1.061405429
49
50 /* ---------------------------------------------------------------------- */
51
PairLJLongCoulLong(LAMMPS * lmp)52 PairLJLongCoulLong::PairLJLongCoulLong(LAMMPS *lmp) : Pair(lmp)
53 {
54 dispersionflag = ewaldflag = pppmflag = 1;
55 respa_enable = 1;
56 writedata = 1;
57 ftable = nullptr;
58 fdisptable = nullptr;
59 qdist = 0.0;
60 cut_respa = nullptr;
61 }
62
63 /* ----------------------------------------------------------------------
64 global settings
65 ------------------------------------------------------------------------- */
66
options(char ** arg,int order)67 void PairLJLongCoulLong::options(char **arg, int order)
68 {
69 const char *option[] = {"long", "cut", "off", nullptr};
70 int i;
71
72 if (!*arg) error->all(FLERR,"Illegal pair_style lj/long/coul/long command");
73 for (i=0; option[i]&&strcmp(arg[0], option[i]); ++i);
74 switch (i) {
75 case 0: ewald_order |= 1<<order; break;
76 case 2: ewald_off |= 1<<order; break;
77 case 1: break;
78 default: error->all(FLERR,"Illegal pair_style lj/long/coul/long command");
79 }
80 }
81
settings(int narg,char ** arg)82 void PairLJLongCoulLong::settings(int narg, char **arg)
83 {
84 if (narg != 3 && narg != 4) error->all(FLERR,"Illegal pair_style command");
85
86 ewald_order = 0;
87 ewald_off = 0;
88
89 options(arg,6);
90 options(++arg,1);
91
92 if (!comm->me && ewald_order == ((1<<1) | (1<<6)))
93 error->warning(FLERR,"Using largest cutoff for lj/long/coul/long");
94 if (!*(++arg))
95 error->all(FLERR,"Cutoffs missing in pair_style lj/long/coul/long");
96 if (!((ewald_order^ewald_off) & (1<<6)))
97 dispersionflag = 0;
98 if (!((ewald_order^ewald_off) & (1<<1)))
99 error->all(FLERR,"Coulomb cut not supported in pair_style lj/long/coul/long");
100 cut_lj_global = utils::numeric(FLERR,*(arg++),false,lmp);
101 if (narg == 4 && ((ewald_order & 0x42) == 0x42))
102 error->all(FLERR,"Only one cutoff allowed when requesting all long");
103 if (narg == 4) cut_coul = utils::numeric(FLERR,*arg,false,lmp);
104 else cut_coul = cut_lj_global;
105
106 if (allocated) {
107 int i,j;
108 for (i = 1; i <= atom->ntypes; i++)
109 for (j = i; j <= atom->ntypes; j++)
110 if (setflag[i][j]) cut_lj[i][j] = cut_lj_global;
111 }
112 }
113
114 /* ----------------------------------------------------------------------
115 free all arrays
116 ------------------------------------------------------------------------- */
117
~PairLJLongCoulLong()118 PairLJLongCoulLong::~PairLJLongCoulLong()
119 {
120 if (allocated) {
121 memory->destroy(setflag);
122 memory->destroy(cutsq);
123
124 memory->destroy(cut_lj_read);
125 memory->destroy(cut_lj);
126 memory->destroy(cut_ljsq);
127 memory->destroy(epsilon_read);
128 memory->destroy(epsilon);
129 memory->destroy(sigma_read);
130 memory->destroy(sigma);
131 memory->destroy(lj1);
132 memory->destroy(lj2);
133 memory->destroy(lj3);
134 memory->destroy(lj4);
135 memory->destroy(offset);
136 }
137 if (ftable) free_tables();
138 if (fdisptable) free_disp_tables();
139 }
140
141 /* ----------------------------------------------------------------------
142 allocate all arrays
143 ------------------------------------------------------------------------- */
144
allocate()145 void PairLJLongCoulLong::allocate()
146 {
147 allocated = 1;
148 int n = atom->ntypes;
149
150 memory->create(setflag,n+1,n+1,"pair:setflag");
151 for (int i = 1; i <= n; i++)
152 for (int j = i; j <= n; j++)
153 setflag[i][j] = 0;
154
155 memory->create(cutsq,n+1,n+1,"pair:cutsq");
156
157 memory->create(cut_lj_read,n+1,n+1,"pair:cut_lj_read");
158 memory->create(cut_lj,n+1,n+1,"pair:cut_lj");
159 memory->create(cut_ljsq,n+1,n+1,"pair:cut_ljsq");
160 memory->create(epsilon_read,n+1,n+1,"pair:epsilon_read");
161 memory->create(epsilon,n+1,n+1,"pair:epsilon");
162 memory->create(sigma_read,n+1,n+1,"pair:sigma_read");
163 memory->create(sigma,n+1,n+1,"pair:sigma");
164 memory->create(lj1,n+1,n+1,"pair:lj1");
165 memory->create(lj2,n+1,n+1,"pair:lj2");
166 memory->create(lj3,n+1,n+1,"pair:lj3");
167 memory->create(lj4,n+1,n+1,"pair:lj4");
168 memory->create(offset,n+1,n+1,"pair:offset");
169 }
170
171 /* ----------------------------------------------------------------------
172 extract protected data from object
173 ------------------------------------------------------------------------- */
174
extract(const char * id,int & dim)175 void *PairLJLongCoulLong::extract(const char *id, int &dim)
176 {
177 const char *ids[] = {
178 "B", "sigma", "epsilon", "ewald_order", "ewald_cut", "ewald_mix",
179 "cut_coul", "cut_LJ", nullptr};
180 void *ptrs[] = {
181 lj4, sigma, epsilon, &ewald_order, &cut_coul, &mix_flag,
182 &cut_coul, &cut_lj_global, nullptr};
183 int i;
184
185 for (i=0; ids[i]&&strcmp(ids[i], id); ++i);
186 if (i <= 2) dim = 2;
187 else dim = 0;
188 return ptrs[i];
189 }
190
191 /* ----------------------------------------------------------------------
192 set coeffs for one or more type pairs
193 ------------------------------------------------------------------------- */
194
coeff(int narg,char ** arg)195 void PairLJLongCoulLong::coeff(int narg, char **arg)
196 {
197 if (narg < 4 || narg > 5) error->all(FLERR,"Incorrect args for pair coefficients");
198 if (!allocated) allocate();
199
200 int ilo,ihi,jlo,jhi;
201 utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error);
202 utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error);
203
204 double epsilon_one = utils::numeric(FLERR,arg[2],false,lmp);
205 double sigma_one = utils::numeric(FLERR,arg[3],false,lmp);
206
207 double cut_lj_one = cut_lj_global;
208 if (narg == 5) cut_lj_one = utils::numeric(FLERR,arg[4],false,lmp);
209
210 int count = 0;
211 for (int i = ilo; i <= ihi; i++) {
212 for (int j = MAX(jlo,i); j <= jhi; j++) {
213 epsilon_read[i][j] = epsilon_one;
214 sigma_read[i][j] = sigma_one;
215 cut_lj_read[i][j] = cut_lj_one;
216 setflag[i][j] = 1;
217 count++;
218 }
219 }
220
221 if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients");
222 }
223
224 /* ----------------------------------------------------------------------
225 init specific to this pair style
226 ------------------------------------------------------------------------- */
227
init_style()228 void PairLJLongCoulLong::init_style()
229 {
230 // require an atom style with charge defined
231
232 if (!atom->q_flag && (ewald_order&(1<<1)))
233 error->all(FLERR,
234 "Invoking coulombic in pair style lj/long/coul/long requires atom attribute q");
235
236 // ensure use of KSpace long-range solver, set two g_ewalds
237
238 if (force->kspace == nullptr)
239 error->all(FLERR,"Pair style requires a KSpace style");
240 if (ewald_order&(1<<1)) g_ewald = force->kspace->g_ewald;
241 if (ewald_order&(1<<6)) g_ewald_6 = force->kspace->g_ewald_6;
242
243 // set rRESPA cutoffs
244
245 if (utils::strmatch(update->integrate_style,"^respa") &&
246 ((Respa *) update->integrate)->level_inner >= 0)
247 cut_respa = ((Respa *) update->integrate)->cutoff;
248 else cut_respa = nullptr;
249
250 // setup force tables
251
252 if (ncoultablebits && (ewald_order&(1<<1))) init_tables(cut_coul,cut_respa);
253 if (ndisptablebits && (ewald_order&(1<<6))) init_tables_disp(cut_lj_global);
254
255 // request regular or rRESPA neighbor lists if neighrequest_flag != 0
256
257 if (force->kspace->neighrequest_flag) {
258 int irequest;
259 int respa = 0;
260
261 if (update->whichflag == 1 && utils::strmatch(update->integrate_style,"^respa")) {
262 if (((Respa *) update->integrate)->level_inner >= 0) respa = 1;
263 if (((Respa *) update->integrate)->level_middle >= 0) respa = 2;
264 }
265
266 irequest = neighbor->request(this,instance_me);
267
268 if (respa >= 1) {
269 neighbor->requests[irequest]->respaouter = 1;
270 neighbor->requests[irequest]->respainner = 1;
271 }
272 if (respa == 2) neighbor->requests[irequest]->respamiddle = 1;
273 }
274
275 cut_coulsq = cut_coul * cut_coul;
276 }
277
278 /* ----------------------------------------------------------------------
279 init for one type pair i,j and corresponding j,i
280 ------------------------------------------------------------------------- */
281
init_one(int i,int j)282 double PairLJLongCoulLong::init_one(int i, int j)
283 {
284 if (setflag[i][j] == 0) {
285 epsilon[i][j] = mix_energy(epsilon_read[i][i],epsilon_read[j][j],
286 sigma_read[i][i],sigma_read[j][j]);
287 sigma[i][j] = mix_distance(sigma_read[i][i],sigma_read[j][j]);
288 if (ewald_order&(1<<6))
289 cut_lj[i][j] = cut_lj_global;
290 else
291 cut_lj[i][j] = mix_distance(cut_lj_read[i][i],cut_lj_read[j][j]);
292 }
293 else {
294 sigma[i][j] = sigma_read[i][j];
295 epsilon[i][j] = epsilon_read[i][j];
296 cut_lj[i][j] = cut_lj_read[i][j];
297 }
298
299 double cut = MAX(cut_lj[i][j], cut_coul + 2.0*qdist);
300 cutsq[i][j] = cut*cut;
301 cut_ljsq[i][j] = cut_lj[i][j] * cut_lj[i][j];
302
303 lj1[i][j] = 48.0 * epsilon[i][j] * pow(sigma[i][j],12.0);
304 lj2[i][j] = 24.0 * epsilon[i][j] * pow(sigma[i][j],6.0);
305 lj3[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],12.0);
306 lj4[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],6.0);
307
308 // check interior rRESPA cutoff
309
310 if (cut_respa && MIN(cut_lj[i][j],cut_coul) < cut_respa[3])
311 error->all(FLERR,"Pair cutoff < Respa interior cutoff");
312
313 if (offset_flag && (cut_lj[i][j] > 0.0)) {
314 double ratio = sigma[i][j] / cut_lj[i][j];
315 offset[i][j] = 4.0 * epsilon[i][j] * (pow(ratio,12.0) - pow(ratio,6.0));
316 } else offset[i][j] = 0.0;
317
318 cutsq[j][i] = cutsq[i][j];
319 cut_ljsq[j][i] = cut_ljsq[i][j];
320 lj1[j][i] = lj1[i][j];
321 lj2[j][i] = lj2[i][j];
322 lj3[j][i] = lj3[i][j];
323 lj4[j][i] = lj4[i][j];
324 offset[j][i] = offset[i][j];
325
326 return cut;
327 }
328
329 /* ----------------------------------------------------------------------
330 proc 0 writes to restart file
331 ------------------------------------------------------------------------- */
332
write_restart(FILE * fp)333 void PairLJLongCoulLong::write_restart(FILE *fp)
334 {
335 write_restart_settings(fp);
336
337 int i,j;
338 for (i = 1; i <= atom->ntypes; i++)
339 for (j = i; j <= atom->ntypes; j++) {
340 fwrite(&setflag[i][j],sizeof(int),1,fp);
341 if (setflag[i][j]) {
342 fwrite(&epsilon_read[i][j],sizeof(double),1,fp);
343 fwrite(&sigma_read[i][j],sizeof(double),1,fp);
344 fwrite(&cut_lj_read[i][j],sizeof(double),1,fp);
345 }
346 }
347 }
348
349 /* ----------------------------------------------------------------------
350 proc 0 reads from restart file, bcasts
351 ------------------------------------------------------------------------- */
352
read_restart(FILE * fp)353 void PairLJLongCoulLong::read_restart(FILE *fp)
354 {
355 read_restart_settings(fp);
356
357 allocate();
358
359 int i,j;
360 int me = comm->me;
361 for (i = 1; i <= atom->ntypes; i++)
362 for (j = i; j <= atom->ntypes; j++) {
363 if (me == 0) utils::sfread(FLERR,&setflag[i][j],sizeof(int),1,fp,nullptr,error);
364 MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world);
365 if (setflag[i][j]) {
366 if (me == 0) {
367 utils::sfread(FLERR,&epsilon_read[i][j],sizeof(double),1,fp,nullptr,error);
368 utils::sfread(FLERR,&sigma_read[i][j],sizeof(double),1,fp,nullptr,error);
369 utils::sfread(FLERR,&cut_lj_read[i][j],sizeof(double),1,fp,nullptr,error);
370 }
371 MPI_Bcast(&epsilon_read[i][j],1,MPI_DOUBLE,0,world);
372 MPI_Bcast(&sigma_read[i][j],1,MPI_DOUBLE,0,world);
373 MPI_Bcast(&cut_lj_read[i][j],1,MPI_DOUBLE,0,world);
374 }
375 }
376 }
377
378 /* ----------------------------------------------------------------------
379 proc 0 writes to restart file
380 ------------------------------------------------------------------------- */
381
write_restart_settings(FILE * fp)382 void PairLJLongCoulLong::write_restart_settings(FILE *fp)
383 {
384 fwrite(&cut_lj_global,sizeof(double),1,fp);
385 fwrite(&cut_coul,sizeof(double),1,fp);
386 fwrite(&offset_flag,sizeof(int),1,fp);
387 fwrite(&mix_flag,sizeof(int),1,fp);
388 fwrite(&ncoultablebits,sizeof(int),1,fp);
389 fwrite(&tabinner,sizeof(double),1,fp);
390 fwrite(&ewald_order,sizeof(int),1,fp);
391 fwrite(&dispersionflag,sizeof(int),1,fp);
392 }
393
394 /* ----------------------------------------------------------------------
395 proc 0 reads from restart file, bcasts
396 ------------------------------------------------------------------------- */
397
read_restart_settings(FILE * fp)398 void PairLJLongCoulLong::read_restart_settings(FILE *fp)
399 {
400 if (comm->me == 0) {
401 utils::sfread(FLERR,&cut_lj_global,sizeof(double),1,fp,nullptr,error);
402 utils::sfread(FLERR,&cut_coul,sizeof(double),1,fp,nullptr,error);
403 utils::sfread(FLERR,&offset_flag,sizeof(int),1,fp,nullptr,error);
404 utils::sfread(FLERR,&mix_flag,sizeof(int),1,fp,nullptr,error);
405 utils::sfread(FLERR,&ncoultablebits,sizeof(int),1,fp,nullptr,error);
406 utils::sfread(FLERR,&tabinner,sizeof(double),1,fp,nullptr,error);
407 utils::sfread(FLERR,&ewald_order,sizeof(int),1,fp,nullptr,error);
408 utils::sfread(FLERR,&dispersionflag,sizeof(int),1,fp,nullptr,error);
409 }
410 MPI_Bcast(&cut_lj_global,1,MPI_DOUBLE,0,world);
411 MPI_Bcast(&cut_coul,1,MPI_DOUBLE,0,world);
412 MPI_Bcast(&offset_flag,1,MPI_INT,0,world);
413 MPI_Bcast(&mix_flag,1,MPI_INT,0,world);
414 MPI_Bcast(&ncoultablebits,1,MPI_INT,0,world);
415 MPI_Bcast(&tabinner,1,MPI_DOUBLE,0,world);
416 MPI_Bcast(&ewald_order,1,MPI_INT,0,world);
417 MPI_Bcast(&dispersionflag,1,MPI_INT,0,world);
418 }
419
420 /* ----------------------------------------------------------------------
421 proc 0 writes to data file
422 ------------------------------------------------------------------------- */
423
write_data(FILE * fp)424 void PairLJLongCoulLong::write_data(FILE *fp)
425 {
426 for (int i = 1; i <= atom->ntypes; i++)
427 fmt::print(fp,"{} {} {}\n",i,epsilon_read[i][i],sigma_read[i][i]);
428 }
429
430 /* ----------------------------------------------------------------------
431 proc 0 writes all pairs to data file. must use the "mixed" parameters.
432 also must not write out cutoff for lj = long
433 ------------------------------------------------------------------------- */
434
write_data_all(FILE * fp)435 void PairLJLongCoulLong::write_data_all(FILE *fp)
436 {
437 for (int i = 1; i <= atom->ntypes; i++) {
438 for (int j = i; j <= atom->ntypes; j++) {
439 if (ewald_order & (1<<6)) {
440 fmt::print(fp,"{} {} {} {}\n",i,j,
441 epsilon[i][j],sigma[i][j]);
442 } else {
443 fmt::print(fp,"{} {} {} {} {}\n",i,j,
444 epsilon[i][j],sigma[i][j],cut_lj[i][j]);
445 }
446 }
447 }
448 }
449
450 /* ----------------------------------------------------------------------
451 compute pair interactions
452 ------------------------------------------------------------------------- */
453
compute(int eflag,int vflag)454 void PairLJLongCoulLong::compute(int eflag, int vflag)
455 {
456 double evdwl,ecoul,fpair;
457 evdwl = ecoul = 0.0;
458 ev_init(eflag,vflag);
459
460 double **x = atom->x, *x0 = x[0];
461 double **f = atom->f, *f0 = f[0], *fi = f0;
462 double *q = atom->q;
463 int *type = atom->type;
464 int nlocal = atom->nlocal;
465 double *special_coul = force->special_coul;
466 double *special_lj = force->special_lj;
467 int newton_pair = force->newton_pair;
468 double qqrd2e = force->qqrd2e;
469
470 int i, j, order1 = ewald_order&(1<<1), order6 = ewald_order&(1<<6);
471 int *ineigh, *ineighn, *jneigh, *jneighn, typei, typej, ni;
472 double qi = 0.0, qri = 0.0;
473 double *cutsqi, *cut_ljsqi, *lj1i, *lj2i, *lj3i, *lj4i, *offseti;
474 double rsq, r2inv, force_coul, force_lj;
475 double g2 = g_ewald_6*g_ewald_6, g6 = g2*g2*g2, g8 = g6*g2;
476 double xi[3], d[3];
477
478 ineighn = (ineigh = list->ilist)+list->inum;
479
480 for (; ineigh<ineighn; ++ineigh) { // loop over my atoms
481 i = *ineigh; fi = f0+3*i;
482 if (order1) qri = (qi = q[i])*qqrd2e; // initialize constants
483 offseti = offset[typei = type[i]];
484 lj1i = lj1[typei]; lj2i = lj2[typei]; lj3i = lj3[typei]; lj4i = lj4[typei];
485 cutsqi = cutsq[typei]; cut_ljsqi = cut_ljsq[typei];
486 memcpy(xi, x0+(i+(i<<1)), 3*sizeof(double));
487 jneighn = (jneigh = list->firstneigh[i])+list->numneigh[i];
488
489 for (; jneigh<jneighn; ++jneigh) { // loop over neighbors
490 j = *jneigh;
491 ni = sbmask(j);
492 j &= NEIGHMASK;
493
494 { double *xj = x0+(j+(j<<1));
495 d[0] = xi[0] - xj[0]; // pair vector
496 d[1] = xi[1] - xj[1];
497 d[2] = xi[2] - xj[2]; }
498
499 if ((rsq = dot3(d, d)) >= cutsqi[typej = type[j]]) continue;
500 r2inv = 1.0/rsq;
501
502 if (order1 && (rsq < cut_coulsq)) { // coulombic
503 if (!ncoultablebits || rsq <= tabinnersq) { // series real space
504 double r = sqrt(rsq), x = g_ewald*r;
505 double s = qri*q[j], t = 1.0/(1.0+EWALD_P*x);
506 if (ni == 0) {
507 s *= g_ewald*exp(-x*x);
508 force_coul = (t *= ((((t*A5+A4)*t+A3)*t+A2)*t+A1)*s/x)+EWALD_F*s;
509 if (eflag) ecoul = t;
510 } else { // special case
511 r = s*(1.0-special_coul[ni])/r; s *= g_ewald*exp(-x*x);
512 force_coul = (t *= ((((t*A5+A4)*t+A3)*t+A2)*t+A1)*s/x)+EWALD_F*s-r;
513 if (eflag) ecoul = t-r;
514 }
515 } else { // table real space
516 union_int_float_t t;
517 t.f = rsq;
518 const int k = (t.i & ncoulmask)>>ncoulshiftbits;
519 double f = (rsq-rtable[k])*drtable[k], qiqj = qi*q[j];
520 if (ni == 0) {
521 force_coul = qiqj*(ftable[k]+f*dftable[k]);
522 if (eflag) ecoul = qiqj*(etable[k]+f*detable[k]);
523 } else { // special case
524 t.f = (1.0-special_coul[ni])*(ctable[k]+f*dctable[k]);
525 force_coul = qiqj*(ftable[k]+f*dftable[k]-t.f);
526 if (eflag) ecoul = qiqj*(etable[k]+f*detable[k]-t.f);
527 }
528 }
529 } else force_coul = ecoul = 0.0;
530
531 if (rsq < cut_ljsqi[typej]) { // lj
532 if (order6) { // long-range lj
533 if (!ndisptablebits || rsq <= tabinnerdispsq) { // series real space
534 double rn = r2inv*r2inv*r2inv;
535 double x2 = g2*rsq, a2 = 1.0/x2;
536 x2 = a2*exp(-x2)*lj4i[typej];
537 if (ni == 0) {
538 force_lj =
539 (rn*=rn)*lj1i[typej]-g8*(((6.0*a2+6.0)*a2+3.0)*a2+1.0)*x2*rsq;
540 if (eflag)
541 evdwl = rn*lj3i[typej]-g6*((a2+1.0)*a2+0.5)*x2;
542 } else { // special case
543 double f = special_lj[ni], t = rn*(1.0-f);
544 force_lj = f*(rn *= rn)*lj1i[typej]-
545 g8*(((6.0*a2+6.0)*a2+3.0)*a2+1.0)*x2*rsq+t*lj2i[typej];
546 if (eflag)
547 evdwl = f*rn*lj3i[typej]-g6*((a2+1.0)*a2+0.5)*x2+t*lj4i[typej];
548 }
549 } else { // table real space
550 union_int_float_t disp_t;
551 disp_t.f = rsq;
552 const int disp_k = (disp_t.i & ndispmask)>>ndispshiftbits;
553 double f_disp = (rsq-rdisptable[disp_k])*drdisptable[disp_k];
554 double rn = r2inv*r2inv*r2inv;
555 if (ni == 0) {
556 force_lj = (rn*=rn)*lj1i[typej]-(fdisptable[disp_k]+f_disp*dfdisptable[disp_k])*lj4i[typej];
557 if (eflag) evdwl = rn*lj3i[typej]-(edisptable[disp_k]+f_disp*dedisptable[disp_k])*lj4i[typej];
558 } else { // special case
559 double f = special_lj[ni], t = rn*(1.0-f);
560 force_lj = f*(rn *= rn)*lj1i[typej]-(fdisptable[disp_k]+f_disp*dfdisptable[disp_k])*lj4i[typej]+t*lj2i[typej];
561 if (eflag) evdwl = f*rn*lj3i[typej]-(edisptable[disp_k]+f_disp*dedisptable[disp_k])*lj4i[typej]+t*lj4i[typej];
562 }
563 }
564 }
565 else { // cut lj
566 double rn = r2inv*r2inv*r2inv;
567 if (ni == 0) {
568 force_lj = rn*(rn*lj1i[typej]-lj2i[typej]);
569 if (eflag) evdwl = rn*(rn*lj3i[typej]-lj4i[typej])-offseti[typej];
570 }
571 else { // special case
572 double f = special_lj[ni];
573 force_lj = f*rn*(rn*lj1i[typej]-lj2i[typej]);
574 if (eflag)
575 evdwl = f * (rn*(rn*lj3i[typej]-lj4i[typej])-offseti[typej]);
576 }
577 }
578 }
579
580 else force_lj = evdwl = 0.0;
581
582 fpair = (force_coul+force_lj)*r2inv;
583
584 if (newton_pair || j < nlocal) {
585 double *fj = f0+(j+(j<<1)), f;
586 fi[0] += f = d[0]*fpair; fj[0] -= f;
587 fi[1] += f = d[1]*fpair; fj[1] -= f;
588 fi[2] += f = d[2]*fpair; fj[2] -= f;
589 }
590 else {
591 fi[0] += d[0]*fpair;
592 fi[1] += d[1]*fpair;
593 fi[2] += d[2]*fpair;
594 }
595
596 if (evflag) ev_tally(i,j,nlocal,newton_pair,
597 evdwl,ecoul,fpair,d[0],d[1],d[2]);
598 }
599 }
600
601 if (vflag_fdotr) virial_fdotr_compute();
602 }
603
604 /* ---------------------------------------------------------------------- */
605
compute_inner()606 void PairLJLongCoulLong::compute_inner()
607 {
608 double rsq, r2inv, force_coul = 0.0, force_lj, fpair;
609
610 int *type = atom->type;
611 int nlocal = atom->nlocal;
612 double *x0 = atom->x[0], *f0 = atom->f[0], *fi = f0, *q = atom->q;
613 double *special_coul = force->special_coul;
614 double *special_lj = force->special_lj;
615 int newton_pair = force->newton_pair;
616 double qqrd2e = force->qqrd2e;
617
618
619 double cut_out_on = cut_respa[0];
620 double cut_out_off = cut_respa[1];
621
622
623 double cut_out_diff = cut_out_off - cut_out_on;
624 double cut_out_on_sq = cut_out_on*cut_out_on;
625 double cut_out_off_sq = cut_out_off*cut_out_off;
626
627 int *ineigh, *ineighn, *jneigh, *jneighn, typei, typej, ni;
628 int i, j, order1 = (ewald_order|(ewald_off^-1))&(1<<1);
629 double qri, *cut_ljsqi, *lj1i, *lj2i;
630 double xi[3], d[3];
631
632 ineighn = (ineigh = list->ilist_inner)+list->inum_inner;
633 for (; ineigh<ineighn; ++ineigh) { // loop over my atoms
634 i = *ineigh; fi = f0+3*i;
635 memcpy(xi, x0+(i+(i<<1)), 3*sizeof(double));
636 cut_ljsqi = cut_ljsq[typei = type[i]];
637 lj1i = lj1[typei]; lj2i = lj2[typei];
638 jneighn = (jneigh = list->firstneigh_inner[i])+list->numneigh_inner[i];
639 for (; jneigh<jneighn; ++jneigh) { // loop over neighbors
640 j = *jneigh;
641 ni = sbmask(j);
642 j &= NEIGHMASK;
643
644 { double *xj = x0+(j+(j<<1));
645 d[0] = xi[0] - xj[0]; // pair vector
646 d[1] = xi[1] - xj[1];
647 d[2] = xi[2] - xj[2]; }
648
649 if ((rsq = dot3(d, d)) >= cut_out_off_sq) continue;
650 r2inv = 1.0/rsq;
651
652 if (order1 && (rsq < cut_coulsq)) { // coulombic
653 qri = qqrd2e*q[i];
654 force_coul = ni == 0 ?
655 qri*q[j]*sqrt(r2inv) : qri*q[j]*sqrt(r2inv)*special_coul[ni];
656 }
657
658 if (rsq < cut_ljsqi[typej = type[j]]) { // lennard-jones
659 double rn = r2inv*r2inv*r2inv;
660 force_lj = ni == 0 ?
661 rn*(rn*lj1i[typej]-lj2i[typej]) :
662 rn*(rn*lj1i[typej]-lj2i[typej])*special_lj[ni];
663 }
664 else force_lj = 0.0;
665
666 fpair = (force_coul + force_lj) * r2inv;
667
668 if (rsq > cut_out_on_sq) { // switching
669 double rsw = (sqrt(rsq) - cut_out_on)/cut_out_diff;
670 fpair *= 1.0 + rsw*rsw*(2.0*rsw-3.0);
671 }
672
673 if (newton_pair || j < nlocal) { // force update
674 double *fj = f0+(j+(j<<1)), f;
675 fi[0] += f = d[0]*fpair; fj[0] -= f;
676 fi[1] += f = d[1]*fpair; fj[1] -= f;
677 fi[2] += f = d[2]*fpair; fj[2] -= f;
678 }
679 else {
680 fi[0] += d[0]*fpair;
681 fi[1] += d[1]*fpair;
682 fi[2] += d[2]*fpair;
683 }
684 }
685 }
686 }
687
688 /* ---------------------------------------------------------------------- */
689
compute_middle()690 void PairLJLongCoulLong::compute_middle()
691 {
692 double rsq, r2inv, force_coul = 0.0, force_lj, fpair;
693
694 int *type = atom->type;
695 int nlocal = atom->nlocal;
696 double *x0 = atom->x[0], *f0 = atom->f[0], *fi = f0, *q = atom->q;
697 double *special_coul = force->special_coul;
698 double *special_lj = force->special_lj;
699 int newton_pair = force->newton_pair;
700 double qqrd2e = force->qqrd2e;
701
702 double cut_in_off = cut_respa[0];
703 double cut_in_on = cut_respa[1];
704 double cut_out_on = cut_respa[2];
705 double cut_out_off = cut_respa[3];
706
707 double cut_in_diff = cut_in_on - cut_in_off;
708 double cut_out_diff = cut_out_off - cut_out_on;
709 double cut_in_off_sq = cut_in_off*cut_in_off;
710 double cut_in_on_sq = cut_in_on*cut_in_on;
711 double cut_out_on_sq = cut_out_on*cut_out_on;
712 double cut_out_off_sq = cut_out_off*cut_out_off;
713
714 int *ineigh, *ineighn, *jneigh, *jneighn, typei, typej, ni;
715 int i, j, order1 = (ewald_order|(ewald_off^-1))&(1<<1);
716 double qri, *cut_ljsqi, *lj1i, *lj2i;
717 double xi[3], d[3];
718
719 ineighn = (ineigh = list->ilist_middle)+list->inum_middle;
720
721 for (; ineigh<ineighn; ++ineigh) { // loop over my atoms
722 i = *ineigh; fi = f0+3*i;
723 if (order1) qri = qqrd2e*q[i];
724 memcpy(xi, x0+(i+(i<<1)), 3*sizeof(double));
725 cut_ljsqi = cut_ljsq[typei = type[i]];
726 lj1i = lj1[typei]; lj2i = lj2[typei];
727 jneighn = (jneigh = list->firstneigh_middle[i])+list->numneigh_middle[i];
728
729 for (; jneigh<jneighn; ++jneigh) {
730 j = *jneigh;
731 ni = sbmask(j);
732 j &= NEIGHMASK;
733
734 { double *xj = x0+(j+(j<<1));
735 d[0] = xi[0] - xj[0]; // pair vector
736 d[1] = xi[1] - xj[1];
737 d[2] = xi[2] - xj[2]; }
738
739 if ((rsq = dot3(d, d)) >= cut_out_off_sq) continue;
740 if (rsq <= cut_in_off_sq) continue;
741 r2inv = 1.0/rsq;
742
743 if (order1 && (rsq < cut_coulsq)) // coulombic
744 force_coul = ni == 0 ?
745 qri*q[j]*sqrt(r2inv) : qri*q[j]*sqrt(r2inv)*special_coul[ni];
746
747 if (rsq < cut_ljsqi[typej = type[j]]) { // lennard-jones
748 double rn = r2inv*r2inv*r2inv;
749 force_lj = ni == 0 ?
750 rn*(rn*lj1i[typej]-lj2i[typej]) :
751 rn*(rn*lj1i[typej]-lj2i[typej])*special_lj[ni];
752 }
753 else force_lj = 0.0;
754
755 fpair = (force_coul + force_lj) * r2inv;
756
757 if (rsq < cut_in_on_sq) { // switching
758 double rsw = (sqrt(rsq) - cut_in_off)/cut_in_diff;
759 fpair *= rsw*rsw*(3.0 - 2.0*rsw);
760 }
761 if (rsq > cut_out_on_sq) {
762 double rsw = (sqrt(rsq) - cut_out_on)/cut_out_diff;
763 fpair *= 1.0 + rsw*rsw*(2.0*rsw-3.0);
764 }
765
766 if (newton_pair || j < nlocal) { // force update
767 double *fj = f0+(j+(j<<1)), f;
768 fi[0] += f = d[0]*fpair; fj[0] -= f;
769 fi[1] += f = d[1]*fpair; fj[1] -= f;
770 fi[2] += f = d[2]*fpair; fj[2] -= f;
771 }
772 else {
773 fi[0] += d[0]*fpair;
774 fi[1] += d[1]*fpair;
775 fi[2] += d[2]*fpair;
776 }
777 }
778 }
779 }
780
781 /* ---------------------------------------------------------------------- */
782
compute_outer(int eflag,int vflag)783 void PairLJLongCoulLong::compute_outer(int eflag, int vflag)
784 {
785 double evdwl,ecoul,fvirial,fpair;
786 evdwl = ecoul = 0.0;
787 ev_init(eflag,vflag);
788
789 double **x = atom->x, *x0 = x[0];
790 double **f = atom->f, *f0 = f[0], *fi = f0;
791 double *q = atom->q;
792 int *type = atom->type;
793 int nlocal = atom->nlocal;
794 double *special_coul = force->special_coul;
795 double *special_lj = force->special_lj;
796 int newton_pair = force->newton_pair;
797 double qqrd2e = force->qqrd2e;
798
799 int i, j, order1 = ewald_order&(1<<1), order6 = ewald_order&(1<<6);
800 int *ineigh, *ineighn, *jneigh, *jneighn, typei, typej, ni, respa_flag;
801 double qi = 0.0, qri = 0.0;
802 double *cutsqi, *cut_ljsqi, *lj1i, *lj2i, *lj3i, *lj4i, *offseti;
803 double rsq, r2inv, force_coul, force_lj;
804 double g2 = g_ewald_6*g_ewald_6, g6 = g2*g2*g2, g8 = g6*g2;
805 double respa_lj = 0.0, respa_coul = 0.0, frespa = 0.0;
806 double xi[3], d[3];
807
808 double cut_in_off = cut_respa[2];
809 double cut_in_on = cut_respa[3];
810
811 double cut_in_diff = cut_in_on - cut_in_off;
812 double cut_in_off_sq = cut_in_off*cut_in_off;
813 double cut_in_on_sq = cut_in_on*cut_in_on;
814
815 ineighn = (ineigh = list->ilist)+list->inum;
816
817 for (; ineigh<ineighn; ++ineigh) { // loop over my atoms
818 i = *ineigh; fi = f0+3*i;
819 if (order1) qri = (qi = q[i])*qqrd2e; // initialize constants
820 offseti = offset[typei = type[i]];
821 lj1i = lj1[typei]; lj2i = lj2[typei]; lj3i = lj3[typei]; lj4i = lj4[typei];
822 cutsqi = cutsq[typei]; cut_ljsqi = cut_ljsq[typei];
823 memcpy(xi, x0+(i+(i<<1)), 3*sizeof(double));
824 jneighn = (jneigh = list->firstneigh[i])+list->numneigh[i];
825
826 for (; jneigh<jneighn; ++jneigh) { // loop over neighbors
827 j = *jneigh;
828 ni = sbmask(j);
829 j &= NEIGHMASK;
830
831 { double *xj = x0+(j+(j<<1));
832 d[0] = xi[0] - xj[0]; // pair vector
833 d[1] = xi[1] - xj[1];
834 d[2] = xi[2] - xj[2]; }
835
836 if ((rsq = dot3(d, d)) >= cutsqi[typej = type[j]]) continue;
837 r2inv = 1.0/rsq;
838
839 frespa = 1.0; // check whether and how to compute respa corrections
840 respa_coul = 0;
841 respa_lj = 0;
842 respa_flag = rsq < cut_in_on_sq ? 1 : 0;
843 if (respa_flag && (rsq > cut_in_off_sq)) {
844 double rsw = (sqrt(rsq)-cut_in_off)/cut_in_diff;
845 frespa = 1-rsw*rsw*(3.0-2.0*rsw);
846 }
847
848 if (order1 && (rsq < cut_coulsq)) { // coulombic
849 if (!ncoultablebits || rsq <= tabinnersq) { // series real space
850 double r = sqrt(rsq), s = qri*q[j];
851 if (respa_flag) // correct for respa
852 respa_coul = ni == 0 ? frespa*s/r : frespa*s/r*special_coul[ni];
853 double x = g_ewald*r, t = 1.0/(1.0+EWALD_P*x);
854 if (ni == 0) {
855 s *= g_ewald*exp(-x*x);
856 force_coul = (t *= ((((t*A5+A4)*t+A3)*t+A2)*t+A1)*s/x)+EWALD_F*s-respa_coul;
857 if (eflag) ecoul = t;
858 }
859 else { // correct for special
860 r = s*(1.0-special_coul[ni])/r; s *= g_ewald*exp(-x*x);
861 force_coul = (t *= ((((t*A5+A4)*t+A3)*t+A2)*t+A1)*s/x)+EWALD_F*s-r-respa_coul;
862 if (eflag) ecoul = t-r;
863 }
864 } // table real space
865 else {
866 if (respa_flag) {
867 double r = sqrt(rsq), s = qri*q[j];
868 respa_coul = ni == 0 ? frespa*s/r : frespa*s/r*special_coul[ni];
869 }
870 union_int_float_t t;
871 t.f = rsq;
872 const int k = (t.i & ncoulmask) >> ncoulshiftbits;
873 double f = (rsq-rtable[k])*drtable[k], qiqj = qi*q[j];
874 if (ni == 0) {
875 force_coul = qiqj*(ftable[k]+f*dftable[k]);
876 if (eflag) ecoul = qiqj*(etable[k]+f*detable[k]);
877 }
878 else { // correct for special
879 t.f = (1.0-special_coul[ni])*(ctable[k]+f*dctable[k]);
880 force_coul = qiqj*(ftable[k]+f*dftable[k]-t.f);
881 if (eflag) {
882 t.f = (1.0-special_coul[ni])*(ptable[k]+f*dptable[k]);
883 ecoul = qiqj*(etable[k]+f*detable[k]-t.f);
884 }
885 }
886 }
887 }
888
889 else force_coul = respa_coul = ecoul = 0.0;
890
891 if (rsq < cut_ljsqi[typej]) { // lennard-jones
892 double rn = r2inv*r2inv*r2inv;
893 if (respa_flag) respa_lj = ni == 0 ? // correct for respa
894 frespa*rn*(rn*lj1i[typej]-lj2i[typej]) :
895 frespa*rn*(rn*lj1i[typej]-lj2i[typej])*special_lj[ni];
896 if (order6) { // long-range form
897 if (!ndisptablebits || rsq <= tabinnerdispsq) {
898 double x2 = g2*rsq, a2 = 1.0/x2;
899 x2 = a2*exp(-x2)*lj4i[typej];
900 if (ni == 0) {
901 force_lj =
902 (rn*=rn)*lj1i[typej]-g8*(((6.0*a2+6.0)*a2+3.0)*a2+1.0)*x2*rsq-respa_lj;
903 if (eflag) evdwl = rn*lj3i[typej]-g6*((a2+1.0)*a2+0.5)*x2;
904 }
905 else { // correct for special
906 double f = special_lj[ni], t = rn*(1.0-f);
907 force_lj = f*(rn *= rn)*lj1i[typej]-
908 g8*(((6.0*a2+6.0)*a2+3.0)*a2+1.0)*x2*rsq+t*lj2i[typej]-respa_lj;
909 if (eflag)
910 evdwl = f*rn*lj3i[typej]-g6*((a2+1.0)*a2+0.5)*x2+t*lj4i[typej];
911 }
912 }
913 else { // table real space
914 union_int_float_t disp_t;
915 disp_t.f = rsq;
916 const int disp_k = (disp_t.i & ndispmask)>>ndispshiftbits;
917 double f_disp = (rsq-rdisptable[disp_k])*drdisptable[disp_k];
918 double rn = r2inv*r2inv*r2inv;
919 if (ni == 0) {
920 force_lj = (rn*=rn)*lj1i[typej]-(fdisptable[disp_k]+f_disp*dfdisptable[disp_k])*lj4i[typej]-respa_lj;
921 if (eflag) evdwl = rn*lj3i[typej]-(edisptable[disp_k]+f_disp*dedisptable[disp_k])*lj4i[typej];
922 }
923 else { // special case
924 double f = special_lj[ni], t = rn*(1.0-f);
925 force_lj = f*(rn *= rn)*lj1i[typej]-(fdisptable[disp_k]+f_disp*dfdisptable[disp_k])*lj4i[typej]+t*lj2i[typej]-respa_lj;
926 if (eflag) evdwl = f*rn*lj3i[typej]-(edisptable[disp_k]+f_disp*dedisptable[disp_k])*lj4i[typej]+t*lj4i[typej];
927 }
928 }
929 }
930 else { // cut form
931 if (ni == 0) {
932 force_lj = rn*(rn*lj1i[typej]-lj2i[typej])-respa_lj;
933 if (eflag) evdwl = rn*(rn*lj3i[typej]-lj4i[typej])-offseti[typej];
934 }
935 else { // correct for special
936 double f = special_lj[ni];
937 force_lj = f*rn*(rn*lj1i[typej]-lj2i[typej])-respa_lj;
938 if (eflag)
939 evdwl = f*(rn*(rn*lj3i[typej]-lj4i[typej])-offseti[typej]);
940 }
941 }
942 }
943 else force_lj = respa_lj = evdwl = 0.0;
944
945 fpair = (force_coul+force_lj)*r2inv;
946
947 if (newton_pair || j < nlocal) {
948 double *fj = f0+(j+(j<<1)), f;
949 fi[0] += f = d[0]*fpair; fj[0] -= f;
950 fi[1] += f = d[1]*fpair; fj[1] -= f;
951 fi[2] += f = d[2]*fpair; fj[2] -= f;
952 }
953 else {
954 fi[0] += d[0]*fpair;
955 fi[1] += d[1]*fpair;
956 fi[2] += d[2]*fpair;
957 }
958
959 if (evflag) {
960 fvirial = (force_coul + force_lj + respa_coul + respa_lj)*r2inv;
961 ev_tally(i,j,nlocal,newton_pair,
962 evdwl,ecoul,fvirial,d[0],d[1],d[2]);
963 }
964 }
965 }
966 }
967
968 /* ---------------------------------------------------------------------- */
969
single(int i,int j,int itype,int jtype,double rsq,double factor_coul,double factor_lj,double & fforce)970 double PairLJLongCoulLong::single(int i, int j, int itype, int jtype,
971 double rsq, double factor_coul, double factor_lj,
972 double &fforce)
973 {
974 double r2inv, r6inv, force_coul, force_lj;
975 double g2 = g_ewald_6*g_ewald_6, g6 = g2*g2*g2, g8 = g6*g2, *q = atom->q;
976
977 double eng = 0.0;
978
979 r2inv = 1.0/rsq;
980 if ((ewald_order&2) && (rsq < cut_coulsq)) { // coulombic
981 if (!ncoultablebits || rsq <= tabinnersq) { // series real space
982 double r = sqrt(rsq), x = g_ewald*r;
983 double s = force->qqrd2e*q[i]*q[j], t = 1.0/(1.0+EWALD_P*x);
984 r = s*(1.0-factor_coul)/r; s *= g_ewald*exp(-x*x);
985 force_coul = (t *= ((((t*A5+A4)*t+A3)*t+A2)*t+A1)*s/x)+EWALD_F*s-r;
986 eng += t-r;
987 } else { // table real space
988 union_int_float_t t;
989 t.f = rsq;
990 const int k = (t.i & ncoulmask) >> ncoulshiftbits;
991 double f = (rsq-rtable[k])*drtable[k], qiqj = q[i]*q[j];
992 t.f = (1.0-factor_coul)*(ctable[k]+f*dctable[k]);
993 force_coul = qiqj*(ftable[k]+f*dftable[k]-t.f);
994 eng += qiqj*(etable[k]+f*detable[k]-t.f);
995 }
996 } else force_coul = 0.0;
997
998 if (rsq < cut_ljsq[itype][jtype]) { // lennard-jones
999 r6inv = r2inv*r2inv*r2inv;
1000 if (ewald_order&64) { // long-range
1001 double x2 = g2*rsq, a2 = 1.0/x2, t = r6inv*(1.0-factor_lj);
1002 x2 = a2*exp(-x2)*lj4[itype][jtype];
1003 force_lj = factor_lj*(r6inv *= r6inv)*lj1[itype][jtype]-
1004 g8*(((6.0*a2+6.0)*a2+3.0)*a2+1.0)*x2*rsq+t*lj2[itype][jtype];
1005 eng += factor_lj*r6inv*lj3[itype][jtype]-
1006 g6*((a2+1.0)*a2+0.5)*x2+t*lj4[itype][jtype];
1007 } else { // cut
1008 force_lj = factor_lj*r6inv*(lj1[itype][jtype]*r6inv-lj2[itype][jtype]);
1009 eng += factor_lj*(r6inv*(r6inv*lj3[itype][jtype]-
1010 lj4[itype][jtype])-offset[itype][jtype]);
1011 }
1012 } else force_lj = 0.0;
1013
1014 fforce = (force_coul+force_lj)*r2inv;
1015 return eng;
1016 }
1017