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