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: Eduardo Bringa (LLNL)
17 ------------------------------------------------------------------------- */
18 
19 #include "pair_buck_coul_cut.h"
20 
21 #include "atom.h"
22 #include "comm.h"
23 #include "error.h"
24 #include "force.h"
25 #include "math_const.h"
26 #include "memory.h"
27 #include "neigh_list.h"
28 #include "neighbor.h"
29 
30 #include <cmath>
31 #include <cstring>
32 
33 using namespace LAMMPS_NS;
34 using namespace MathConst;
35 
36 /* ---------------------------------------------------------------------- */
37 
PairBuckCoulCut(LAMMPS * lmp)38 PairBuckCoulCut::PairBuckCoulCut(LAMMPS *lmp) : Pair(lmp)
39 {
40   writedata = 1;
41 }
42 
43 /* ---------------------------------------------------------------------- */
44 
~PairBuckCoulCut()45 PairBuckCoulCut::~PairBuckCoulCut()
46 {
47   if (copymode) return;
48 
49   if (allocated) {
50     memory->destroy(setflag);
51     memory->destroy(cutsq);
52 
53     memory->destroy(cut_lj);
54     memory->destroy(cut_ljsq);
55     memory->destroy(cut_coul);
56     memory->destroy(cut_coulsq);
57     memory->destroy(a);
58     memory->destroy(rho);
59     memory->destroy(c);
60     memory->destroy(rhoinv);
61     memory->destroy(buck1);
62     memory->destroy(buck2);
63     memory->destroy(offset);
64   }
65 }
66 
67 /* ---------------------------------------------------------------------- */
68 
compute(int eflag,int vflag)69 void PairBuckCoulCut::compute(int eflag, int vflag)
70 {
71   int i,j,ii,jj,inum,jnum,itype,jtype;
72   double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,evdwl,ecoul,fpair;
73   double rsq,r2inv,r6inv,r,rexp,forcecoul,forcebuck,factor_coul,factor_lj;
74   int *ilist,*jlist,*numneigh,**firstneigh;
75 
76   evdwl = ecoul = 0.0;
77   ev_init(eflag,vflag);
78 
79   double **x = atom->x;
80   double **f = atom->f;
81   double *q = atom->q;
82   int *type = atom->type;
83   int nlocal = atom->nlocal;
84   double *special_coul = force->special_coul;
85   double *special_lj = force->special_lj;
86   int newton_pair = force->newton_pair;
87   double qqrd2e = force->qqrd2e;
88 
89   inum = list->inum;
90   ilist = list->ilist;
91   numneigh = list->numneigh;
92   firstneigh = list->firstneigh;
93 
94   // loop over neighbors of my atoms
95 
96   for (ii = 0; ii < inum; ii++) {
97     i = ilist[ii];
98     qtmp = q[i];
99     xtmp = x[i][0];
100     ytmp = x[i][1];
101     ztmp = x[i][2];
102     itype = type[i];
103     jlist = firstneigh[i];
104     jnum = numneigh[i];
105 
106     for (jj = 0; jj < jnum; jj++) {
107       j = jlist[jj];
108       factor_lj = special_lj[sbmask(j)];
109       factor_coul = special_coul[sbmask(j)];
110       j &= NEIGHMASK;
111 
112       delx = xtmp - x[j][0];
113       dely = ytmp - x[j][1];
114       delz = ztmp - x[j][2];
115       rsq = delx*delx + dely*dely + delz*delz;
116       jtype = type[j];
117 
118       if (rsq < cutsq[itype][jtype]) {
119         r2inv = 1.0/rsq;
120         r = sqrt(rsq);
121 
122         if (rsq < cut_coulsq[itype][jtype])
123           forcecoul = qqrd2e * qtmp*q[j]/r;
124         else forcecoul = 0.0;
125 
126         if (rsq < cut_ljsq[itype][jtype]) {
127           r6inv = r2inv*r2inv*r2inv;
128           rexp = exp(-r*rhoinv[itype][jtype]);
129           forcebuck = buck1[itype][jtype]*r*rexp - buck2[itype][jtype]*r6inv;
130         } else forcebuck = 0.0;
131 
132         fpair = (factor_coul*forcecoul + factor_lj*forcebuck) * r2inv;
133 
134         f[i][0] += delx*fpair;
135         f[i][1] += dely*fpair;
136         f[i][2] += delz*fpair;
137         if (newton_pair || j < nlocal) {
138           f[j][0] -= delx*fpair;
139           f[j][1] -= dely*fpair;
140           f[j][2] -= delz*fpair;
141         }
142 
143         if (eflag) {
144           if (rsq < cut_coulsq[itype][jtype])
145             ecoul = factor_coul * qqrd2e * qtmp*q[j]/r;
146           else ecoul = 0.0;
147           if (rsq < cut_ljsq[itype][jtype]) {
148             evdwl = a[itype][jtype]*rexp - c[itype][jtype]*r6inv -
149               offset[itype][jtype];
150             evdwl *= factor_lj;
151           } else evdwl = 0.0;
152         }
153 
154         if (evflag) ev_tally(i,j,nlocal,newton_pair,
155                              evdwl,ecoul,fpair,delx,dely,delz);
156       }
157     }
158   }
159 
160   if (vflag_fdotr) virial_fdotr_compute();
161 }
162 
163 /* ----------------------------------------------------------------------
164    allocate all arrays
165 ------------------------------------------------------------------------- */
166 
allocate()167 void PairBuckCoulCut::allocate()
168 {
169   allocated = 1;
170   int n = atom->ntypes;
171 
172   memory->create(setflag,n+1,n+1,"pair:setflag");
173 
174   for (int i = 1; i <= n; i++)
175     for (int j = i; j <= n; j++)
176       setflag[i][j] = 0;
177 
178   memory->create(cutsq,n+1,n+1,"pair:cutsq");
179 
180   memory->create(cut_lj,n+1,n+1,"pair:cut_lj");
181   memory->create(cut_ljsq,n+1,n+1,"pair:cut_ljsq");
182   memory->create(cut_coul,n+1,n+1,"pair:cut_coul");
183   memory->create(cut_coulsq,n+1,n+1,"pair:cut_coulsq");
184   memory->create(a,n+1,n+1,"pair:a");
185   memory->create(rho,n+1,n+1,"pair:rho");
186   memory->create(c,n+1,n+1,"pair:c");
187   memory->create(rhoinv,n+1,n+1,"pair:rhoinv");
188   memory->create(buck1,n+1,n+1,"pair:buck1");
189   memory->create(buck2,n+1,n+1,"pair:buck2");
190   memory->create(offset,n+1,n+1,"pair:offset");
191 }
192 
193 /* ----------------------------------------------------------------------
194    global settings
195 ------------------------------------------------------------------------- */
196 
settings(int narg,char ** arg)197 void PairBuckCoulCut::settings(int narg, char **arg)
198 {
199   if (narg < 1 || narg > 2) error->all(FLERR,"Illegal pair_style command");
200 
201   cut_lj_global = utils::numeric(FLERR,arg[0],false,lmp);
202   if (narg == 1) cut_coul_global = cut_lj_global;
203   else cut_coul_global = utils::numeric(FLERR,arg[1],false,lmp);
204 
205   // reset cutoffs that have been explicitly set
206 
207   if (allocated) {
208     int i,j;
209     for (i = 1; i <= atom->ntypes; i++)
210       for (j = i; j <= atom->ntypes; j++)
211         if (setflag[i][j]) {
212           cut_lj[i][j] = cut_lj_global;
213           cut_coul[i][j] = cut_coul_global;
214         }
215   }
216 }
217 
218 /* ----------------------------------------------------------------------
219    set coeffs for one or more type pairs
220 ------------------------------------------------------------------------- */
221 
coeff(int narg,char ** arg)222 void PairBuckCoulCut::coeff(int narg, char **arg)
223 {
224   if (narg < 5 || narg > 7)
225     error->all(FLERR,"Incorrect args for pair coefficients");
226   if (!allocated) allocate();
227 
228   int ilo,ihi,jlo,jhi;
229   utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error);
230   utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error);
231 
232   double a_one = utils::numeric(FLERR,arg[2],false,lmp);
233   double rho_one = utils::numeric(FLERR,arg[3],false,lmp);
234   if (rho_one <= 0) error->all(FLERR,"Incorrect args for pair coefficients");
235   double c_one = utils::numeric(FLERR,arg[4],false,lmp);
236 
237   double cut_lj_one = cut_lj_global;
238   double cut_coul_one = cut_coul_global;
239   if (narg >= 6) cut_coul_one = cut_lj_one = utils::numeric(FLERR,arg[5],false,lmp);
240   if (narg == 7) cut_coul_one = utils::numeric(FLERR,arg[6],false,lmp);
241 
242   int count = 0;
243   for (int i = ilo; i <= ihi; i++) {
244     for (int j = MAX(jlo,i); j <= jhi; j++) {
245       a[i][j] = a_one;
246       rho[i][j] = rho_one;
247       c[i][j] = c_one;
248       cut_lj[i][j] = cut_lj_one;
249       cut_coul[i][j] = cut_coul_one;
250       setflag[i][j] = 1;
251       count++;
252     }
253   }
254 
255   if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients");
256 }
257 
258 /* ----------------------------------------------------------------------
259    init specific to this pair style
260 ------------------------------------------------------------------------- */
261 
init_style()262 void PairBuckCoulCut::init_style()
263 {
264   if (!atom->q_flag)
265     error->all(FLERR,"Pair style buck/coul/cut requires atom attribute q");
266 
267   neighbor->request(this,instance_me);
268 }
269 
270 /* ----------------------------------------------------------------------
271    init for one type pair i,j and corresponding j,i
272 ------------------------------------------------------------------------- */
273 
init_one(int i,int j)274 double PairBuckCoulCut::init_one(int i, int j)
275 {
276   if (setflag[i][j] == 0) error->all(FLERR,"All pair coeffs are not set");
277 
278   double cut = MAX(cut_lj[i][j],cut_coul[i][j]);
279   cut_ljsq[i][j] = cut_lj[i][j] * cut_lj[i][j];
280   cut_coulsq[i][j] = cut_coul[i][j] * cut_coul[i][j];
281 
282   rhoinv[i][j] = 1.0/rho[i][j];
283   buck1[i][j] = a[i][j]/rho[i][j];
284   buck2[i][j] = 6.0*c[i][j];
285 
286   if (offset_flag && (cut_lj[i][j] > 0.0)) {
287     double rexp = exp(-cut_lj[i][j]/rho[i][j]);
288     offset[i][j] = a[i][j]*rexp - c[i][j]/pow(cut_lj[i][j],6.0);
289   } else offset[i][j] = 0.0;
290 
291   cut_ljsq[j][i] = cut_ljsq[i][j];
292   cut_coulsq[j][i] = cut_coulsq[i][j];
293   a[j][i] = a[i][j];
294   c[j][i] = c[i][j];
295   rhoinv[j][i] = rhoinv[i][j];
296   buck1[j][i] = buck1[i][j];
297   buck2[j][i] = buck2[i][j];
298   offset[j][i] = offset[i][j];
299 
300   // compute I,J contribution to long-range tail correction
301   // count total # of atoms of type I and J via Allreduce
302 
303   if (tail_flag) {
304     int *type = atom->type;
305     int nlocal = atom->nlocal;
306 
307     double count[2],all[2];
308     count[0] = count[1] = 0.0;
309     for (int k = 0; k < nlocal; k++) {
310       if (type[k] == i) count[0] += 1.0;
311       if (type[k] == j) count[1] += 1.0;
312     }
313     MPI_Allreduce(count,all,2,MPI_DOUBLE,MPI_SUM,world);
314 
315     double rho1 = rho[i][j];
316     double rho2 = rho1*rho1;
317     double rho3 = rho2*rho1;
318     double rc = cut_lj[i][j];
319     double rc2 = rc*rc;
320     double rc3 = rc2*rc;
321     etail_ij = 2.0*MY_PI*all[0]*all[1]*
322       (a[i][j]*exp(-rc/rho1)*rho1*(rc2 + 2.0*rho1*rc + 2.0*rho2) -
323        c[i][j]/(3.0*rc3));
324     ptail_ij = (-1/3.0)*2.0*MY_PI*all[0]*all[1]*
325       (-a[i][j]*exp(-rc/rho1)*
326        (rc3 + 3.0*rho1*rc2 + 6.0*rho2*rc + 6.0*rho3) + 2.0*c[i][j]/rc3);
327   }
328 
329   return cut;
330 }
331 
332 /* ----------------------------------------------------------------------
333   proc 0 writes to restart file
334 ------------------------------------------------------------------------- */
335 
write_restart(FILE * fp)336 void PairBuckCoulCut::write_restart(FILE *fp)
337 {
338   write_restart_settings(fp);
339 
340   int i,j;
341   for (i = 1; i <= atom->ntypes; i++)
342     for (j = i; j <= atom->ntypes; j++) {
343       fwrite(&setflag[i][j],sizeof(int),1,fp);
344       if (setflag[i][j]) {
345         fwrite(&a[i][j],sizeof(double),1,fp);
346         fwrite(&rho[i][j],sizeof(double),1,fp);
347         fwrite(&c[i][j],sizeof(double),1,fp);
348         fwrite(&cut_lj[i][j],sizeof(double),1,fp);
349         fwrite(&cut_coul[i][j],sizeof(double),1,fp);
350       }
351     }
352 }
353 
354 /* ----------------------------------------------------------------------
355   proc 0 reads from restart file, bcasts
356 ------------------------------------------------------------------------- */
357 
read_restart(FILE * fp)358 void PairBuckCoulCut::read_restart(FILE *fp)
359 {
360   read_restart_settings(fp);
361 
362   allocate();
363 
364   int i,j;
365   int me = comm->me;
366   for (i = 1; i <= atom->ntypes; i++)
367     for (j = i; j <= atom->ntypes; j++) {
368       if (me == 0) utils::sfread(FLERR,&setflag[i][j],sizeof(int),1,fp,nullptr,error);
369       MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world);
370       if (setflag[i][j]) {
371         if (me == 0) {
372           utils::sfread(FLERR,&a[i][j],sizeof(double),1,fp,nullptr,error);
373           utils::sfread(FLERR,&rho[i][j],sizeof(double),1,fp,nullptr,error);
374           utils::sfread(FLERR,&c[i][j],sizeof(double),1,fp,nullptr,error);
375           utils::sfread(FLERR,&cut_lj[i][j],sizeof(double),1,fp,nullptr,error);
376           utils::sfread(FLERR,&cut_coul[i][j],sizeof(double),1,fp,nullptr,error);
377         }
378         MPI_Bcast(&a[i][j],1,MPI_DOUBLE,0,world);
379         MPI_Bcast(&rho[i][j],1,MPI_DOUBLE,0,world);
380         MPI_Bcast(&c[i][j],1,MPI_DOUBLE,0,world);
381         MPI_Bcast(&cut_lj[i][j],1,MPI_DOUBLE,0,world);
382         MPI_Bcast(&cut_coul[i][j],1,MPI_DOUBLE,0,world);
383       }
384     }
385 }
386 
387 /* ----------------------------------------------------------------------
388   proc 0 writes to restart file
389 ------------------------------------------------------------------------- */
390 
write_restart_settings(FILE * fp)391 void PairBuckCoulCut::write_restart_settings(FILE *fp)
392 {
393   fwrite(&cut_lj_global,sizeof(double),1,fp);
394   fwrite(&cut_coul_global,sizeof(double),1,fp);
395   fwrite(&offset_flag,sizeof(int),1,fp);
396   fwrite(&mix_flag,sizeof(int),1,fp);
397   fwrite(&tail_flag,sizeof(int),1,fp);
398 }
399 
400 /* ----------------------------------------------------------------------
401   proc 0 reads from restart file, bcasts
402 ------------------------------------------------------------------------- */
403 
read_restart_settings(FILE * fp)404 void PairBuckCoulCut::read_restart_settings(FILE *fp)
405 {
406   if (comm->me == 0) {
407     utils::sfread(FLERR,&cut_lj_global,sizeof(double),1,fp,nullptr,error);
408     utils::sfread(FLERR,&cut_coul_global,sizeof(double),1,fp,nullptr,error);
409     utils::sfread(FLERR,&offset_flag,sizeof(int),1,fp,nullptr,error);
410     utils::sfread(FLERR,&mix_flag,sizeof(int),1,fp,nullptr,error);
411     utils::sfread(FLERR,&tail_flag,sizeof(int),1,fp,nullptr,error);
412   }
413   MPI_Bcast(&cut_lj_global,1,MPI_DOUBLE,0,world);
414   MPI_Bcast(&cut_coul_global,1,MPI_DOUBLE,0,world);
415   MPI_Bcast(&offset_flag,1,MPI_INT,0,world);
416   MPI_Bcast(&mix_flag,1,MPI_INT,0,world);
417   MPI_Bcast(&tail_flag,1,MPI_INT,0,world);
418 }
419 
420 /* ----------------------------------------------------------------------
421    proc 0 writes to data file
422 ------------------------------------------------------------------------- */
423 
write_data(FILE * fp)424 void PairBuckCoulCut::write_data(FILE *fp)
425 {
426   for (int i = 1; i <= atom->ntypes; i++)
427     fprintf(fp,"%d %g %g %g\n",i,a[i][i],rho[i][i],c[i][i]);
428 }
429 
430 /* ----------------------------------------------------------------------
431    proc 0 writes all pairs to data file
432 ------------------------------------------------------------------------- */
433 
write_data_all(FILE * fp)434 void PairBuckCoulCut::write_data_all(FILE *fp)
435 {
436   for (int i = 1; i <= atom->ntypes; i++)
437     for (int j = i; j <= atom->ntypes; j++)
438       fprintf(fp,"%d %d %g %g %g %g %g\n",i,j,
439               a[i][j],rho[i][j],c[i][j],cut_lj[i][j],cut_coul[i][j]);
440 }
441 
442 /* ---------------------------------------------------------------------- */
443 
single(int i,int j,int itype,int jtype,double rsq,double factor_coul,double factor_lj,double & fforce)444 double PairBuckCoulCut::single(int i, int j, int itype, int jtype,
445                                double rsq,
446                                double factor_coul, double factor_lj,
447                                double &fforce)
448 {
449   double r2inv,r6inv,r,rexp,forcecoul,forcebuck,phicoul,phibuck;
450 
451   r2inv = 1.0/rsq;
452   if (rsq < cut_coulsq[itype][jtype])
453     forcecoul = force->qqrd2e * atom->q[i]*atom->q[j]*sqrt(r2inv);
454   else forcecoul = 0.0;
455   if (rsq < cut_ljsq[itype][jtype]) {
456     r6inv = r2inv*r2inv*r2inv;
457     r = sqrt(rsq);
458     rexp = exp(-r*rhoinv[itype][jtype]);
459     forcebuck = buck1[itype][jtype]*r*rexp - buck2[itype][jtype]*r6inv;
460   } else forcebuck = 0.0;
461   fforce = (factor_coul*forcecoul + factor_lj*forcebuck) * r2inv;
462 
463   double eng = 0.0;
464   if (rsq < cut_coulsq[itype][jtype]) {
465     phicoul = force->qqrd2e * atom->q[i]*atom->q[j]*sqrt(r2inv);
466     eng += factor_coul*phicoul;
467   }
468   if (rsq < cut_ljsq[itype][jtype]) {
469     phibuck = a[itype][jtype]*rexp - c[itype][jtype]*r6inv -
470       offset[itype][jtype];
471     eng += factor_lj*phibuck;
472   }
473   return eng;
474 }
475 
476 /* ---------------------------------------------------------------------- */
477 
extract(const char * str,int & dim)478 void *PairBuckCoulCut::extract(const char *str, int &dim)
479 {
480   dim = 2;
481   if (strcmp(str,"a") == 0) return (void *) a;
482   if (strcmp(str,"c") == 0) return (void *) c;
483   return nullptr;
484 }
485