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:  Evangelos Voyiatzis (Royal DSM)
17  * ------------------------------------------------------------------------- */
18 
19 #include "pair_coul_slater_long.h"
20 
21 #include "atom.h"
22 #include "comm.h"
23 #include "error.h"
24 #include "force.h"
25 #include "kspace.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 
35 #define EWALD_F   1.12837917
36 #define EWALD_P   0.3275911
37 #define A1        0.254829592
38 #define A2       -0.284496736
39 #define A3        1.421413741
40 #define A4       -1.453152027
41 #define A5        1.061405429
42 
43 /* ---------------------------------------------------------------------- */
44 
PairCoulSlaterLong(LAMMPS * lmp)45 PairCoulSlaterLong::PairCoulSlaterLong(LAMMPS *lmp) : Pair(lmp)
46 {
47   ewaldflag = pppmflag = 1;
48   //ftable = nullptr;
49   qdist = 0.0;
50 }
51 
52 /* ---------------------------------------------------------------------- */
53 
~PairCoulSlaterLong()54 PairCoulSlaterLong::~PairCoulSlaterLong()
55 {
56   if (!copymode) {
57     if (allocated) {
58       memory->destroy(setflag);
59       memory->destroy(cutsq);
60 
61       memory->destroy(scale);
62     }
63   }
64 }
65 
66 /* ---------------------------------------------------------------------- */
67 
compute(int eflag,int vflag)68 void PairCoulSlaterLong::compute(int eflag, int vflag)
69 {
70   int i,j,ii,jj,inum,jnum,itype,jtype;
71   double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,ecoul,fpair;
72   double r,r2inv,forcecoul,factor_coul;
73   double grij,expm2,prefactor,t,erfc;
74   int *ilist,*jlist,*numneigh,**firstneigh;
75   double rsq;
76   double slater_term;
77 
78   ev_init(eflag,vflag);
79   ecoul = 0.0;
80 
81   double **x = atom->x;
82   double **f = atom->f;
83   double *q = atom->q;
84   int *type = atom->type;
85   int nlocal = atom->nlocal;
86   double *special_coul = force->special_coul;
87   int newton_pair = force->newton_pair;
88   double qqrd2e = force->qqrd2e;
89 
90   inum = list->inum;
91   ilist = list->ilist;
92   numneigh = list->numneigh;
93   firstneigh = list->firstneigh;
94 
95   // loop over neighbors of my atoms
96 
97   for (ii = 0; ii < inum; ii++) {
98     i = ilist[ii];
99     qtmp = q[i];
100     xtmp = x[i][0];
101     ytmp = x[i][1];
102     ztmp = x[i][2];
103     itype = type[i];
104     jlist = firstneigh[i];
105     jnum = numneigh[i];
106 
107     for (jj = 0; jj < jnum; jj++) {
108       j = jlist[jj];
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 < cut_coulsq) {
119         r2inv = 1.0/rsq;
120 //        if (!ncoultablebits || rsq <= tabinnersq) {
121           r = sqrt(rsq);
122           grij = g_ewald * r;
123           expm2 = exp(-grij*grij);
124           t = 1.0 / (1.0 + EWALD_P*grij);
125           erfc = t * (A1+t*(A2+t*(A3+t*(A4+t*A5)))) * expm2;
126           slater_term = exp(-2*r/lamda)*(1 + (2*r/lamda*(1+r/lamda)));
127           prefactor = qqrd2e * scale[itype][jtype] * qtmp*q[j]/r;
128           forcecoul = prefactor * (erfc + EWALD_F*grij*expm2 - slater_term);
129           if (factor_coul < 1.0) forcecoul -= (1.0-factor_coul)*prefactor;
130 /*
131         } else {
132           union_int_float_t rsq_lookup;
133           rsq_lookup.f = rsq;
134           itable = rsq_lookup.i & ncoulmask;
135           itable >>= ncoulshiftbits;
136           fraction = (rsq_lookup.f - rtable[itable]) * drtable[itable];
137           table = ftable[itable] + fraction*dftable[itable];
138           forcecoul = scale[itype][jtype] * qtmp*q[j] * table;
139           if (factor_coul < 1.0) {
140             table = ctable[itable] + fraction*dctable[itable];
141             prefactor = scale[itype][jtype] * qtmp*q[j] * table;
142             forcecoul -= (1.0-factor_coul)*prefactor;
143           }
144         }
145 */
146 
147         fpair = forcecoul * r2inv;
148 
149         f[i][0] += delx*fpair;
150         f[i][1] += dely*fpair;
151         f[i][2] += delz*fpair;
152         if (newton_pair || j < nlocal) {
153           f[j][0] -= delx*fpair;
154           f[j][1] -= dely*fpair;
155           f[j][2] -= delz*fpair;
156         }
157 
158         if (eflag) {
159 //          if (!ncoultablebits || rsq <= tabinnersq)
160             ecoul = prefactor*(erfc - (1 + r/lamda)*exp(-2*r/lamda));
161 /*
162           else {
163             table = etable[itable] + fraction*detable[itable];
164             ecoul = scale[itype][jtype] * qtmp*q[j] * table;
165           }
166 */
167           if (factor_coul < 1.0) ecoul -= (1.0-factor_coul)*prefactor;
168         }
169 
170         if (evflag) ev_tally(i,j,nlocal,newton_pair,
171                              0.0,ecoul,fpair,delx,dely,delz);
172       }
173     }
174   }
175 
176   if (vflag_fdotr) virial_fdotr_compute();
177 }
178 
179 /* ----------------------------------------------------------------------
180    allocate all arrays
181 ------------------------------------------------------------------------- */
182 
allocate()183 void PairCoulSlaterLong::allocate()
184 {
185   allocated = 1;
186   int n = atom->ntypes;
187 
188   memory->create(setflag,n+1,n+1,"pair:setflag");
189   for (int i = 1; i <= n; i++)
190     for (int j = i; j <= n; j++)
191       setflag[i][j] = 0;
192 
193   memory->create(cutsq,n+1,n+1,"pair:cutsq");
194 
195   memory->create(scale,n+1,n+1,"pair:scale");
196 }
197 
198 /* ----------------------------------------------------------------------
199    global settings
200 ------------------------------------------------------------------------- */
201 
settings(int narg,char ** arg)202 void PairCoulSlaterLong::settings(int narg, char **arg)
203 {
204   if (narg != 2) error->all(FLERR,"Illegal pair_style command");
205 
206   lamda = utils::numeric(FLERR,arg[0],false,lmp);
207   cut_coul = utils::numeric(FLERR,arg[1],false,lmp);
208 }
209 
210 /* ----------------------------------------------------------------------
211    set coeffs for one or more type pairs
212 ------------------------------------------------------------------------- */
213 
coeff(int narg,char ** arg)214 void PairCoulSlaterLong::coeff(int narg, char **arg)
215 {
216   if (narg != 2) error->all(FLERR,"Incorrect args for pair coefficients");
217   if (!allocated) allocate();
218 
219   int ilo,ihi,jlo,jhi;
220   utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error);
221   utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error);
222 
223   int count = 0;
224   for (int i = ilo; i <= ihi; i++) {
225     for (int j = MAX(jlo,i); j <= jhi; j++) {
226       scale[i][j] = 1.0;
227       setflag[i][j] = 1;
228       count++;
229     }
230   }
231 
232   if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients");
233 }
234 
235 /* ----------------------------------------------------------------------
236    init specific to this pair style
237 ------------------------------------------------------------------------- */
238 
init_style()239 void PairCoulSlaterLong::init_style()
240 {
241   if (!atom->q_flag)
242     error->all(FLERR,"Pair style coul/slater/long requires atom attribute q");
243 
244   neighbor->request(this,instance_me);
245 
246   cut_coulsq = cut_coul * cut_coul;
247 
248   // insure use of KSpace long-range solver, set g_ewald
249 
250  if (force->kspace == nullptr)
251     error->all(FLERR,"Pair style requires a KSpace style");
252   g_ewald = force->kspace->g_ewald;
253 
254   // setup force tables
255 
256   // if (ncoultablebits) init_tables(cut_coul,nullptr);
257 }
258 
259 /* ----------------------------------------------------------------------
260    init for one type pair i,j and corresponding j,i
261 ------------------------------------------------------------------------- */
262 
init_one(int i,int j)263 double PairCoulSlaterLong::init_one(int i, int j)
264 {
265   scale[j][i] = scale[i][j];
266   return cut_coul+2.0*qdist;
267 }
268 
269 /* ----------------------------------------------------------------------
270   proc 0 writes to restart file
271 ------------------------------------------------------------------------- */
272 
write_restart(FILE * fp)273 void PairCoulSlaterLong::write_restart(FILE *fp)
274 {
275   write_restart_settings(fp);
276 
277   for (int i = 1; i <= atom->ntypes; i++)
278     for (int j = i; j <= atom->ntypes; j++) {
279       fwrite(&setflag[i][j],sizeof(int),1,fp);
280       if (setflag[i][j])
281         fwrite(&scale[i][j],sizeof(double),1,fp);
282     }
283 }
284 
285 /* ----------------------------------------------------------------------
286   proc 0 reads from restart file, bcasts
287 ------------------------------------------------------------------------- */
288 
read_restart(FILE * fp)289 void PairCoulSlaterLong::read_restart(FILE *fp)
290 {
291   read_restart_settings(fp);
292 
293   allocate();
294 
295   int i,j;
296   int me = comm->me;
297   for (i = 1; i <= atom->ntypes; i++)
298     for (j = i; j <= atom->ntypes; j++) {
299       if (me == 0) fread(&setflag[i][j],sizeof(int),1,fp);
300       MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world);
301       if (setflag[i][j]) {
302         if (me == 0) fread(&scale[i][j],sizeof(double),1,fp);
303         MPI_Bcast(&scale[i][j],1,MPI_DOUBLE,0,world);
304       }
305     }
306 }
307 
308 /* ----------------------------------------------------------------------
309   proc 0 writes to restart file
310 ------------------------------------------------------------------------- */
311 
write_restart_settings(FILE * fp)312 void PairCoulSlaterLong::write_restart_settings(FILE *fp)
313 {
314   fwrite(&cut_coul,sizeof(double),1,fp);
315   fwrite(&lamda,sizeof(double),1,fp);
316   fwrite(&offset_flag,sizeof(int),1,fp);
317   fwrite(&mix_flag,sizeof(int),1,fp);
318   //fwrite(&ncoultablebits,sizeof(int),1,fp);
319   //fwrite(&tabinner,sizeof(double),1,fp);
320 }
321 
322 /* ----------------------------------------------------------------------
323   proc 0 reads from restart file, bcasts
324 ------------------------------------------------------------------------- */
325 
read_restart_settings(FILE * fp)326 void PairCoulSlaterLong::read_restart_settings(FILE *fp)
327 {
328   if (comm->me == 0) {
329     fread(&cut_coul,sizeof(double),1,fp);
330     fread(&lamda,sizeof(double),1,fp);
331     fread(&offset_flag,sizeof(int),1,fp);
332     fread(&mix_flag,sizeof(int),1,fp);
333     //fread(&ncoultablebits,sizeof(int),1,fp);
334     //fread(&tabinner,sizeof(double),1,fp);
335   }
336   MPI_Bcast(&cut_coul,1,MPI_DOUBLE,0,world);
337   MPI_Bcast(&lamda,1,MPI_DOUBLE,0,world);
338   MPI_Bcast(&offset_flag,1,MPI_INT,0,world);
339   MPI_Bcast(&mix_flag,1,MPI_INT,0,world);
340   //MPI_Bcast(&ncoultablebits,1,MPI_INT,0,world);
341   //MPI_Bcast(&tabinner,1,MPI_DOUBLE,0,world);
342 }
343 
344 /* ---------------------------------------------------------------------- */
345 
single(int i,int j,int,int,double rsq,double factor_coul,double,double & fforce)346 double PairCoulSlaterLong::single(int i, int j, int /*itype*/, int /*jtype*/,
347                             double rsq,
348                             double factor_coul, double /*factor_lj*/,
349                             double &fforce)
350 {
351   double r2inv,r,grij,expm2,t,erfc,prefactor;
352   double slater_term;
353 //  double fraction,table;
354   double forcecoul,phicoul;
355 //  int itable;
356 
357   r2inv = 1.0/rsq;
358 //  if (!ncoultablebits || rsq <= tabinnersq) {
359     r = sqrt(rsq);
360     grij = g_ewald * r;
361     expm2 = exp(-grij*grij);
362     t = 1.0 / (1.0 + EWALD_P*grij);
363     erfc = t * (A1+t*(A2+t*(A3+t*(A4+t*A5)))) * expm2;
364     slater_term = exp(-2*r/lamda)*(1 + (2*r/lamda*(1+r/lamda)));
365     prefactor = force->qqrd2e * atom->q[i]*atom->q[j]/r;
366     forcecoul = prefactor * (erfc + EWALD_F*grij*expm2 - slater_term);
367     if (factor_coul < 1.0) forcecoul -= (1.0-factor_coul)*prefactor;
368 /*
369   } else {
370     union_int_float_t rsq_lookup;
371     rsq_lookup.f = rsq;
372     itable = rsq_lookup.i & ncoulmask;
373     itable >>= ncoulshiftbits;
374     fraction = (rsq_lookup.f - rtable[itable]) * drtable[itable];
375     table = ftable[itable] + fraction*dftable[itable];
376     forcecoul = atom->q[i]*atom->q[j] * table;
377     if (factor_coul < 1.0) {
378       table = ctable[itable] + fraction*dctable[itable];
379       prefactor = atom->q[i]*atom->q[j] * table;
380       forcecoul -= (1.0-factor_coul)*prefactor;
381     }
382   }
383 */
384   fforce = forcecoul * r2inv;
385 
386 //  if (!ncoultablebits || rsq <= tabinnersq)
387     phicoul = prefactor*(erfc - (1 + r/lamda)*exp(-2*r/lamda));
388 /*
389   else {
390     table = etable[itable] + fraction*detable[itable];
391     phicoul = atom->q[i]*atom->q[j] * table;
392   }
393 */
394   if (factor_coul < 1.0) phicoul -= (1.0-factor_coul)*prefactor;
395 
396   return phicoul;
397 }
398 
399 /* ---------------------------------------------------------------------- */
400 
extract(const char * str,int & dim)401 void *PairCoulSlaterLong::extract(const char *str, int &dim)
402 {
403   if (strcmp(str,"cut_coul") == 0) {
404     dim = 0;
405     return (void *) &cut_coul;
406   }
407   if (strcmp(str,"lamda") == 0) {
408     dim = 0;
409     return (void *) &lamda;
410   }
411   if (strcmp(str,"scale") == 0) {
412     dim = 2;
413     return (void *) scale;
414   }
415   return nullptr;
416 }
417