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