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: Pavel Elkind (Gothenburg University)
17 ------------------------------------------------------------------------- */
18 
19 #include "pair_lj_cut_tip4p_cut.h"
20 
21 #include <cmath>
22 #include <cstring>
23 #include "atom.h"
24 #include "force.h"
25 #include "neighbor.h"
26 #include "neigh_list.h"
27 #include "domain.h"
28 #include "angle.h"
29 #include "bond.h"
30 #include "comm.h"
31 #include "math_const.h"
32 #include "memory.h"
33 #include "error.h"
34 
35 
36 using namespace LAMMPS_NS;
37 using namespace MathConst;
38 
39 /* ---------------------------------------------------------------------- */
40 
PairLJCutTIP4PCut(LAMMPS * lmp)41 PairLJCutTIP4PCut::PairLJCutTIP4PCut(LAMMPS *lmp) : Pair(lmp)
42 {
43   single_enable = 0;
44   writedata = 1;
45 
46   nmax = 0;
47   hneigh = nullptr;
48   newsite = nullptr;
49 
50   // TIP4P cannot compute virial as F dot r
51   // due to finding bonded H atoms which are not near O atom
52 
53   no_virial_fdotr_compute = 1;
54 }
55 
56 /* ---------------------------------------------------------------------- */
57 
~PairLJCutTIP4PCut()58 PairLJCutTIP4PCut::~PairLJCutTIP4PCut()
59 {
60   if (allocated) {
61     memory->destroy(setflag);
62     memory->destroy(cutsq);
63 
64     memory->destroy(cut_lj);
65     memory->destroy(cut_ljsq);
66     memory->destroy(epsilon);
67     memory->destroy(sigma);
68     memory->destroy(lj1);
69     memory->destroy(lj2);
70     memory->destroy(lj3);
71     memory->destroy(lj4);
72     memory->destroy(offset);
73   }
74 
75   memory->destroy(hneigh);
76   memory->destroy(newsite);
77 }
78 
79 /* ---------------------------------------------------------------------- */
80 
compute(int eflag,int vflag)81 void PairLJCutTIP4PCut::compute(int eflag, int vflag)
82 {
83   int i,j,ii,jj,inum,jnum,itype,jtype;
84   double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,evdwl,ecoul;
85   double rsq,r2inv,r6inv,forcecoul,forcelj,factor_lj,factor_coul;
86   int *ilist,*jlist,*numneigh,**firstneigh;
87 
88   int key;
89   int n,vlist[6];
90   int iH1,iH2,jH1,jH2;
91   double cforce;
92   double fO[3],fH[3],fd[3],v[6];
93   double *x1,*x2,*xH1,*xH2;
94 
95   evdwl = ecoul = 0.0;
96   ev_init(eflag,vflag);
97 
98   // reallocate hneigh & newsite if necessary
99   // initialize hneigh[0] to -1 on steps when reneighboring occurred
100   // initialize hneigh[2] to 0 every step
101 
102   int nlocal = atom->nlocal;
103   int nall = nlocal + atom->nghost;
104 
105   if (atom->nmax > nmax) {
106     nmax = atom->nmax;
107     memory->destroy(hneigh);
108     memory->create(hneigh,nmax,3,"pair:hneigh");
109     memory->destroy(newsite);
110     memory->create(newsite,nmax,3,"pair:newsite");
111   }
112   if (neighbor->ago == 0)
113     for (i = 0; i < nall; i++) hneigh[i][0] = -1;
114   for (i = 0; i < nall; i++) hneigh[i][2] = 0;
115 
116   double **f = atom->f;
117   double **x = atom->x;
118   double *q = atom->q;
119   tagint *tag = atom->tag;
120   int *type = atom->type;
121   double *special_lj = force->special_lj;
122   double *special_coul = force->special_coul;
123   int newton_pair = force->newton_pair;
124   double qqrd2e = force->qqrd2e;
125 
126   inum = list->inum;
127   ilist = list->ilist;
128   numneigh = list->numneigh;
129   firstneigh = list->firstneigh;
130 
131   // loop over neighbors of my atoms
132 
133   for (ii = 0; ii < inum; ii++) {
134     i = ilist[ii];
135     qtmp = q[i];
136     xtmp = x[i][0];
137     ytmp = x[i][1];
138     ztmp = x[i][2];
139     itype = type[i];
140 
141     if (itype == typeO) {
142       if (hneigh[i][0] < 0) {
143         iH1 = atom->map(tag[i] + 1);
144         iH2 = atom->map(tag[i] + 2);
145         if (iH1 == -1 || iH2 == -1)
146           error->one(FLERR,"TIP4P hydrogen is missing");
147         if (atom->type[iH1] != typeH || atom->type[iH2] != typeH)
148           error->one(FLERR,"TIP4P hydrogen has incorrect atom type");
149         // set iH1,iH2 to index of closest image to O
150         iH1 = domain->closest_image(i,iH1);
151         iH2 = domain->closest_image(i,iH2);
152         compute_newsite(x[i],x[iH1],x[iH2],newsite[i]);
153         hneigh[i][0] = iH1;
154         hneigh[i][1] = iH2;
155         hneigh[i][2] = 1;
156 
157       } else {
158         iH1 = hneigh[i][0];
159         iH2 = hneigh[i][1];
160         if (hneigh[i][2] == 0) {
161           hneigh[i][2] = 1;
162           compute_newsite(x[i],x[iH1],x[iH2],newsite[i]);
163         }
164       }
165       x1 = newsite[i];
166     } else x1 = x[i];
167 
168     jlist = firstneigh[i];
169     jnum = numneigh[i];
170 
171     for (jj = 0; jj < jnum; jj++) {
172       j = jlist[jj];
173       factor_lj = special_lj[sbmask(j)];
174       factor_coul = special_coul[sbmask(j)];
175       j &= NEIGHMASK;
176 
177       delx = xtmp - x[j][0];
178       dely = ytmp - x[j][1];
179       delz = ztmp - x[j][2];
180       rsq = delx*delx + dely*dely + delz*delz;
181       jtype = type[j];
182 
183       // LJ interaction based on true rsq
184 
185       if (rsq < cut_ljsq[itype][jtype]) {
186         r2inv = 1.0/rsq;
187         r6inv = r2inv*r2inv*r2inv;
188         forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]);
189         forcelj *= factor_lj * r2inv;
190 
191         f[i][0] += delx*forcelj;
192         f[i][1] += dely*forcelj;
193         f[i][2] += delz*forcelj;
194         f[j][0] -= delx*forcelj;
195         f[j][1] -= dely*forcelj;
196         f[j][2] -= delz*forcelj;
197 
198         if (eflag) {
199           evdwl = r6inv*(lj3[itype][jtype]*r6inv-lj4[itype][jtype]) -
200             offset[itype][jtype];
201           evdwl *= factor_lj;
202         } else evdwl = 0.0;
203 
204         if (evflag) ev_tally(i,j,nlocal,newton_pair,
205                              evdwl,0.0,forcelj,delx,dely,delz);
206       }
207 
208       // adjust rsq and delxyz for off-site O charge(s) if necessary
209       // but only if they are within reach
210 
211       if (rsq < cut_coulsqplus) {
212         if (itype == typeO || jtype == typeO) {
213 
214           // if atom J = water O, set x2 = offset charge site
215           // else x2 = x of atom J
216 
217           if (jtype == typeO) {
218             if (hneigh[j][0] < 0) {
219               jH1 = atom->map(tag[j] + 1);
220               jH2 = atom->map(tag[j] + 2);
221               if (jH1 == -1 || jH2 == -1)
222                 error->one(FLERR,"TIP4P hydrogen is missing");
223               if (atom->type[jH1] != typeH || atom->type[jH2] != typeH)
224                 error->one(FLERR,"TIP4P hydrogen has incorrect atom type");
225               // set jH1,jH2 to closest image to O
226               jH1 = domain->closest_image(j,jH1);
227               jH2 = domain->closest_image(j,jH2);
228               compute_newsite(x[j],x[jH1],x[jH2],newsite[j]);
229               hneigh[j][0] = jH1;
230               hneigh[j][1] = jH2;
231               hneigh[j][2] = 1;
232 
233             } else {
234               jH1 = hneigh[j][0];
235               jH2 = hneigh[j][1];
236               if (hneigh[j][2] == 0) {
237                 hneigh[j][2] = 1;
238                 compute_newsite(x[j],x[jH1],x[jH2],newsite[j]);
239               }
240             }
241             x2 = newsite[j];
242           } else x2 = x[j];
243 
244           delx = x1[0] - x2[0];
245           dely = x1[1] - x2[1];
246           delz = x1[2] - x2[2];
247           rsq = delx*delx + dely*dely + delz*delz;
248         }
249 
250         // Coulombic interaction based on modified rsq
251 
252         if (rsq < cut_coulsq) {
253           r2inv = 1.0 / rsq;
254           forcecoul = qqrd2e * qtmp * q[j] * sqrt(r2inv);
255           cforce = factor_coul * forcecoul * r2inv;
256 
257         // if i,j are not O atoms, force is applied directly;
258         // if i or j are O atoms, force is on fictitious atom & partitioned
259         // force partitioning due to Feenstra, J Comp Chem, 20, 786 (1999)
260         // f_f = fictitious force, fO = f_f (1 - 2 alpha), fH = alpha f_f
261         // preserves total force and torque on water molecule
262         // virial = sum(r x F) where each water's atoms are near xi and xj
263         // vlist stores 2,4,6 atoms whose forces contribute to virial
264 
265           n = 0;
266           key = 0;
267 
268           if (itype != typeO) {
269             f[i][0] += delx * cforce;
270             f[i][1] += dely * cforce;
271             f[i][2] += delz * cforce;
272 
273             if (vflag) {
274               v[0] = x[i][0] * delx * cforce;
275               v[1] = x[i][1] * dely * cforce;
276               v[2] = x[i][2] * delz * cforce;
277               v[3] = x[i][0] * dely * cforce;
278               v[4] = x[i][0] * delz * cforce;
279               v[5] = x[i][1] * delz * cforce;
280             }
281             vlist[n++] = i;
282 
283           } else {
284             key++;
285 
286             fd[0] = delx*cforce;
287             fd[1] = dely*cforce;
288             fd[2] = delz*cforce;
289 
290             fO[0] = fd[0]*(1.0 - alpha);
291             fO[1] = fd[1]*(1.0 - alpha);
292             fO[2] = fd[2]*(1.0 - alpha);
293 
294             fH[0] = 0.5 * alpha * fd[0];
295             fH[1] = 0.5 * alpha * fd[1];
296             fH[2] = 0.5 * alpha * fd[2];
297 
298             f[i][0] += fO[0];
299             f[i][1] += fO[1];
300             f[i][2] += fO[2];
301 
302             f[iH1][0] += fH[0];
303             f[iH1][1] += fH[1];
304             f[iH1][2] += fH[2];
305 
306             f[iH2][0] += fH[0];
307             f[iH2][1] += fH[1];
308             f[iH2][2] += fH[2];
309 
310             if (vflag) {
311               xH1 = x[iH1];
312               xH2 = x[iH2];
313               v[0] = x[i][0]*fO[0] + xH1[0]*fH[0] + xH2[0]*fH[0];
314               v[1] = x[i][1]*fO[1] + xH1[1]*fH[1] + xH2[1]*fH[1];
315               v[2] = x[i][2]*fO[2] + xH1[2]*fH[2] + xH2[2]*fH[2];
316               v[3] = x[i][0]*fO[1] + xH1[0]*fH[1] + xH2[0]*fH[1];
317               v[4] = x[i][0]*fO[2] + xH1[0]*fH[2] + xH2[0]*fH[2];
318               v[5] = x[i][1]*fO[2] + xH1[1]*fH[2] + xH2[1]*fH[2];
319             }
320             vlist[n++] = i;
321             vlist[n++] = iH1;
322             vlist[n++] = iH2;
323           }
324 
325           if (jtype != typeO) {
326             f[j][0] -= delx * cforce;
327             f[j][1] -= dely * cforce;
328             f[j][2] -= delz * cforce;
329 
330             if (vflag) {
331               v[0] -= x[j][0] * delx * cforce;
332               v[1] -= x[j][1] * dely * cforce;
333               v[2] -= x[j][2] * delz * cforce;
334               v[3] -= x[j][0] * dely * cforce;
335               v[4] -= x[j][0] * delz * cforce;
336               v[5] -= x[j][1] * delz * cforce;
337             }
338             vlist[n++] = j;
339 
340           } else {
341             key += 2;
342 
343             fd[0] = -delx*cforce;
344             fd[1] = -dely*cforce;
345             fd[2] = -delz*cforce;
346 
347             fO[0] = fd[0]*(1 - alpha);
348             fO[1] = fd[1]*(1 - alpha);
349             fO[2] = fd[2]*(1 - alpha);
350 
351             fH[0] = 0.5 * alpha * fd[0];
352             fH[1] = 0.5 * alpha * fd[1];
353             fH[2] = 0.5 * alpha * fd[2];
354 
355             f[j][0] += fO[0];
356             f[j][1] += fO[1];
357             f[j][2] += fO[2];
358 
359             f[jH1][0] += fH[0];
360             f[jH1][1] += fH[1];
361             f[jH1][2] += fH[2];
362 
363             f[jH2][0] += fH[0];
364             f[jH2][1] += fH[1];
365             f[jH2][2] += fH[2];
366 
367             if (vflag) {
368               xH1 = x[jH1];
369               xH2 = x[jH2];
370               v[0] += x[j][0]*fO[0] + xH1[0]*fH[0] + xH2[0]*fH[0];
371               v[1] += x[j][1]*fO[1] + xH1[1]*fH[1] + xH2[1]*fH[1];
372               v[2] += x[j][2]*fO[2] + xH1[2]*fH[2] + xH2[2]*fH[2];
373               v[3] += x[j][0]*fO[1] + xH1[0]*fH[1] + xH2[0]*fH[1];
374               v[4] += x[j][0]*fO[2] + xH1[0]*fH[2] + xH2[0]*fH[2];
375               v[5] += x[j][1]*fO[2] + xH1[1]*fH[2] + xH2[1]*fH[2];
376             }
377             vlist[n++] = j;
378             vlist[n++] = jH1;
379             vlist[n++] = jH2;
380           }
381 
382           if (eflag) {
383             ecoul = qqrd2e * qtmp * q[j] * sqrt(r2inv);
384             ecoul *= factor_coul;
385           } else ecoul = 0.0;
386 
387           if (evflag) ev_tally_tip4p(key,vlist,v,ecoul,alpha);
388         }
389       }
390     }
391   }
392 }
393 
394 /* ----------------------------------------------------------------------
395    allocate all arrays
396 ------------------------------------------------------------------------- */
397 
allocate()398 void PairLJCutTIP4PCut::allocate()
399 {
400   allocated = 1;
401   int n = atom->ntypes;
402 
403   memory->create(setflag,n+1,n+1,"pair:setflag");
404   for (int i = 1; i <= n; i++)
405     for (int j = i; j <= n; j++)
406       setflag[i][j] = 0;
407 
408   memory->create(cutsq,n+1,n+1,"pair:cutsq");
409 
410   memory->create(cut_lj,n+1,n+1,"pair:cut_lj");
411   memory->create(cut_ljsq,n+1,n+1,"pair:cut_ljsq");
412   memory->create(epsilon,n+1,n+1,"pair:epsilon");
413   memory->create(sigma,n+1,n+1,"pair:sigma");
414   memory->create(lj1,n+1,n+1,"pair:lj1");
415   memory->create(lj2,n+1,n+1,"pair:lj2");
416   memory->create(lj3,n+1,n+1,"pair:lj3");
417   memory->create(lj4,n+1,n+1,"pair:lj4");
418   memory->create(offset,n+1,n+1,"pair:offset");
419 }
420 
421 /* ----------------------------------------------------------------------
422    global settings
423 ------------------------------------------------------------------------- */
424 
settings(int narg,char ** arg)425 void PairLJCutTIP4PCut::settings(int narg, char **arg)
426 {
427   if (narg < 6 || narg > 7) error->all(FLERR,"Illegal pair_style command");
428 
429   typeO = utils::inumeric(FLERR,arg[0],false,lmp);
430   typeH = utils::inumeric(FLERR,arg[1],false,lmp);
431   typeB = utils::inumeric(FLERR,arg[2],false,lmp);
432   typeA = utils::inumeric(FLERR,arg[3],false,lmp);
433   qdist = utils::numeric(FLERR,arg[4],false,lmp);
434 
435   cut_lj_global = utils::numeric(FLERR,arg[5],false,lmp);
436   if (narg == 6) cut_coul = cut_lj_global;
437   else cut_coul = utils::numeric(FLERR,arg[6],false,lmp);
438 
439   cut_coulsq = cut_coul * cut_coul;
440   cut_coulsqplus = (cut_coul + 2.0*qdist) * (cut_coul + 2.0*qdist);
441 
442   if (allocated) {
443     int i,j;
444     for (i = 1; i <= atom->ntypes; i++)
445       for (j = i; j <= atom->ntypes; j++)
446         if (setflag[i][j]) cut_lj[i][j] = cut_lj_global;
447   }
448 }
449 
450 /* ----------------------------------------------------------------------
451    set coeffs for one or more type pairs
452 ------------------------------------------------------------------------- */
453 
coeff(int narg,char ** arg)454 void PairLJCutTIP4PCut::coeff(int narg, char **arg)
455 {
456   if (narg < 4 || narg > 5)
457     error->all(FLERR,"Incorrect args for pair coefficients");
458   if (!allocated) allocate();
459 
460   int ilo,ihi,jlo,jhi;
461   utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error);
462   utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error);
463 
464   double epsilon_one = utils::numeric(FLERR,arg[2],false,lmp);
465   double sigma_one = utils::numeric(FLERR,arg[3],false,lmp);
466 
467   double cut_lj_one = cut_lj_global;
468   if (narg == 5) cut_lj_one = utils::numeric(FLERR,arg[4],false,lmp);
469 
470   int count = 0;
471   for (int i = ilo; i <= ihi; i++) {
472     for (int j = MAX(jlo,i); j <= jhi; j++) {
473       epsilon[i][j] = epsilon_one;
474       sigma[i][j] = sigma_one;
475       cut_lj[i][j] = cut_lj_one;
476       setflag[i][j] = 1;
477       count++;
478     }
479   }
480 
481   if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients");
482 }
483 
484 /* ----------------------------------------------------------------------
485    init specific to this pair style
486 ------------------------------------------------------------------------- */
487 
init_style()488 void PairLJCutTIP4PCut::init_style()
489 {
490   if (atom->tag_enable == 0)
491     error->all(FLERR,"Pair style lj/cut/tip4p/cut requires atom IDs");
492   if (!force->newton_pair)
493     error->all(FLERR,
494                "Pair style lj/cut/tip4p/cut requires newton pair on");
495   if (!atom->q_flag)
496     error->all(FLERR,
497                "Pair style lj/cut/tip4p/cut requires atom attribute q");
498   if (force->bond == nullptr)
499     error->all(FLERR,"Must use a bond style with TIP4P potential");
500   if (force->angle == nullptr)
501     error->all(FLERR,"Must use an angle style with TIP4P potential");
502 
503   neighbor->request(this,instance_me);
504 
505   // set alpha parameter
506 
507   double theta = force->angle->equilibrium_angle(typeA);
508   double blen = force->bond->equilibrium_distance(typeB);
509   alpha = qdist / (cos(0.5*theta) * blen);
510 }
511 
512 /* ----------------------------------------------------------------------
513    init for one type pair i,j and corresponding j,i
514 ------------------------------------------------------------------------- */
515 
init_one(int i,int j)516 double PairLJCutTIP4PCut::init_one(int i, int j)
517 {
518   if (setflag[i][j] == 0) {
519     epsilon[i][j] = mix_energy(epsilon[i][i],epsilon[j][j],
520                                sigma[i][i],sigma[j][j]);
521     sigma[i][j] = mix_distance(sigma[i][i],sigma[j][j]);
522     cut_lj[i][j] = mix_distance(cut_lj[i][i],cut_lj[j][j]);
523   }
524 
525   // include TIP4P qdist in full cutoff, qdist = 0.0 if not TIP4P
526 
527   double cut = MAX(cut_lj[i][j],cut_coul+2.0*qdist);
528   cut_ljsq[i][j] = cut_lj[i][j] * cut_lj[i][j];
529 
530   lj1[i][j] = 48.0 * epsilon[i][j] * pow(sigma[i][j],12.0);
531   lj2[i][j] = 24.0 * epsilon[i][j] * pow(sigma[i][j],6.0);
532   lj3[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],12.0);
533   lj4[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],6.0);
534 
535   if (offset_flag && (cut_lj[i][j] > 0.0)) {
536     double ratio = sigma[i][j] / cut_lj[i][j];
537     offset[i][j] = 4.0 * epsilon[i][j] * (pow(ratio,12.0) - pow(ratio,6.0));
538   } else offset[i][j] = 0.0;
539 
540   cut_ljsq[j][i] = cut_ljsq[i][j];
541   lj1[j][i] = lj1[i][j];
542   lj2[j][i] = lj2[i][j];
543   lj3[j][i] = lj3[i][j];
544   lj4[j][i] = lj4[i][j];
545   offset[j][i] = offset[i][j];
546 
547   // compute I,J contribution to long-range tail correction
548   // count total # of atoms of type I and J via Allreduce
549 
550   if (tail_flag) {
551     int *type = atom->type;
552     int nlocal = atom->nlocal;
553 
554     double count[2],all[2];
555     count[0] = count[1] = 0.0;
556     for (int k = 0; k < nlocal; k++) {
557       if (type[k] == i) count[0] += 1.0;
558       if (type[k] == j) count[1] += 1.0;
559     }
560     MPI_Allreduce(count,all,2,MPI_DOUBLE,MPI_SUM,world);
561 
562     double sig2 = sigma[i][j]*sigma[i][j];
563     double sig6 = sig2*sig2*sig2;
564     double rc3 = cut_lj[i][j]*cut_lj[i][j]*cut_lj[i][j];
565     double rc6 = rc3*rc3;
566     double rc9 = rc3*rc6;
567     etail_ij = 8.0*MY_PI*all[0]*all[1]*epsilon[i][j] *
568       sig6 * (sig6 - 3.0*rc6) / (9.0*rc9);
569     ptail_ij = 16.0*MY_PI*all[0]*all[1]*epsilon[i][j] *
570       sig6 * (2.0*sig6 - 3.0*rc6) / (9.0*rc9);
571   }
572 
573   // check that LJ epsilon = 0.0 for water H
574   // set LJ cutoff to 0.0 for any interaction involving water H
575   // so LJ term isn't calculated in compute()
576 
577   if ((i == typeH && epsilon[i][i] != 0.0) ||
578       (j == typeH && epsilon[j][j] != 0.0))
579     error->all(FLERR,"Water H epsilon must be 0.0 for "
580                "pair style lj/cut/tip4p/cut");
581 
582   if (i == typeH || j == typeH)
583     cut_ljsq[j][i] = cut_ljsq[i][j] = 0.0;
584 
585   return cut;
586 }
587 
588 /* ----------------------------------------------------------------------
589   proc 0 writes to restart file
590 ------------------------------------------------------------------------- */
591 
write_restart(FILE * fp)592 void PairLJCutTIP4PCut::write_restart(FILE *fp)
593 {
594   write_restart_settings(fp);
595 
596   int i,j;
597   for (i = 1; i <= atom->ntypes; i++) {
598     for (j = i; j <= atom->ntypes; j++) {
599       fwrite(&setflag[i][j],sizeof(int),1,fp);
600       if (setflag[i][j]) {
601         fwrite(&epsilon[i][j],sizeof(double),1,fp);
602         fwrite(&sigma[i][j],sizeof(double),1,fp);
603         fwrite(&cut_lj[i][j],sizeof(double),1,fp);
604       }
605     }
606   }
607 }
608 
609 /* ----------------------------------------------------------------------
610   proc 0 reads from restart file, bcasts
611 ------------------------------------------------------------------------- */
612 
read_restart(FILE * fp)613 void PairLJCutTIP4PCut::read_restart(FILE *fp)
614 {
615   read_restart_settings(fp);
616   allocate();
617 
618   int i,j;
619   int me = comm->me;
620   for (i = 1; i <= atom->ntypes; i++) {
621     for (j = i; j <= atom->ntypes; j++) {
622       if (me == 0) utils::sfread(FLERR,&setflag[i][j],sizeof(int),1,fp,nullptr,error);
623       MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world);
624       if (setflag[i][j]) {
625         if (me == 0) {
626           utils::sfread(FLERR,&epsilon[i][j],sizeof(double),1,fp,nullptr,error);
627           utils::sfread(FLERR,&sigma[i][j],sizeof(double),1,fp,nullptr,error);
628           utils::sfread(FLERR,&cut_lj[i][j],sizeof(double),1,fp,nullptr,error);
629         }
630         MPI_Bcast(&epsilon[i][j],1,MPI_DOUBLE,0,world);
631         MPI_Bcast(&sigma[i][j],1,MPI_DOUBLE,0,world);
632         MPI_Bcast(&cut_lj[i][j],1,MPI_DOUBLE,0,world);
633       }
634     }
635   }
636 }
637 
638 /* ----------------------------------------------------------------------
639   proc 0 writes to restart file
640 ------------------------------------------------------------------------- */
641 
write_restart_settings(FILE * fp)642 void PairLJCutTIP4PCut::write_restart_settings(FILE *fp)
643 {
644   fwrite(&typeO,sizeof(int),1,fp);
645   fwrite(&typeH,sizeof(int),1,fp);
646   fwrite(&typeB,sizeof(int),1,fp);
647   fwrite(&typeA,sizeof(int),1,fp);
648   fwrite(&qdist,sizeof(double),1,fp);
649 
650   fwrite(&cut_lj_global,sizeof(double),1,fp);
651   fwrite(&cut_coul,sizeof(double),1,fp);
652   fwrite(&offset_flag,sizeof(int),1,fp);
653   fwrite(&mix_flag,sizeof(int),1,fp);
654   fwrite(&tail_flag,sizeof(int),1,fp);
655 }
656 
657 /* ----------------------------------------------------------------------
658   proc 0 reads from restart file, bcasts
659 ------------------------------------------------------------------------- */
660 
read_restart_settings(FILE * fp)661 void PairLJCutTIP4PCut::read_restart_settings(FILE *fp)
662 {
663   if (comm->me == 0) {
664     utils::sfread(FLERR,&typeO,sizeof(int),1,fp,nullptr,error);
665     utils::sfread(FLERR,&typeH,sizeof(int),1,fp,nullptr,error);
666     utils::sfread(FLERR,&typeB,sizeof(int),1,fp,nullptr,error);
667     utils::sfread(FLERR,&typeA,sizeof(int),1,fp,nullptr,error);
668     utils::sfread(FLERR,&qdist,sizeof(double),1,fp,nullptr,error);
669 
670     utils::sfread(FLERR,&cut_lj_global,sizeof(double),1,fp,nullptr,error);
671     utils::sfread(FLERR,&cut_coul,sizeof(double),1,fp,nullptr,error);
672     utils::sfread(FLERR,&offset_flag,sizeof(int),1,fp,nullptr,error);
673     utils::sfread(FLERR,&mix_flag,sizeof(int),1,fp,nullptr,error);
674     utils::sfread(FLERR,&tail_flag,sizeof(int),1,fp,nullptr,error);
675   }
676 
677   MPI_Bcast(&typeO,1,MPI_INT,0,world);
678   MPI_Bcast(&typeH,1,MPI_INT,0,world);
679   MPI_Bcast(&typeB,1,MPI_INT,0,world);
680   MPI_Bcast(&typeA,1,MPI_INT,0,world);
681   MPI_Bcast(&qdist,1,MPI_DOUBLE,0,world);
682 
683   MPI_Bcast(&cut_lj_global,1,MPI_DOUBLE,0,world);
684   MPI_Bcast(&cut_coul,1,MPI_DOUBLE,0,world);
685   MPI_Bcast(&offset_flag,1,MPI_INT,0,world);
686   MPI_Bcast(&mix_flag,1,MPI_INT,0,world);
687   MPI_Bcast(&tail_flag,1,MPI_INT,0,world);
688 
689   cut_coulsq = cut_coul * cut_coul;
690   cut_coulsqplus = (cut_coul + 2.0*qdist) * (cut_coul + 2.0*qdist);
691 }
692 
693 /* ----------------------------------------------------------------------
694    proc 0 writes to data file
695 ------------------------------------------------------------------------- */
696 
write_data(FILE * fp)697 void PairLJCutTIP4PCut::write_data(FILE *fp)
698 {
699   for (int i = 1; i <= atom->ntypes; i++)
700     fprintf(fp,"%d %g %g\n",i,epsilon[i][i],sigma[i][i]);
701 }
702 
703 /* ----------------------------------------------------------------------
704    proc 0 writes all pairs to data file
705 ------------------------------------------------------------------------- */
706 
write_data_all(FILE * fp)707 void PairLJCutTIP4PCut::write_data_all(FILE *fp)
708 {
709   for (int i = 1; i <= atom->ntypes; i++)
710     for (int j = i; j <= atom->ntypes; j++)
711       fprintf(fp,"%d %d %g %g %g\n",i,j,epsilon[i][j],sigma[i][j],cut_lj[i][j]);
712 }
713 
714 /* ----------------------------------------------------------------------
715   compute position xM of fictitious charge site for O atom and 2 H atoms
716   return it as xM
717 ------------------------------------------------------------------------- */
718 
compute_newsite(double * xO,double * xH1,double * xH2,double * xM)719 void PairLJCutTIP4PCut::compute_newsite(double *xO,  double *xH1,
720                                         double *xH2, double *xM)
721 {
722   double delx1 = xH1[0] - xO[0];
723   double dely1 = xH1[1] - xO[1];
724   double delz1 = xH1[2] - xO[2];
725 
726   double delx2 = xH2[0] - xO[0];
727   double dely2 = xH2[1] - xO[1];
728   double delz2 = xH2[2] - xO[2];
729 
730   xM[0] = xO[0] + alpha * 0.5 * (delx1 + delx2);
731   xM[1] = xO[1] + alpha * 0.5 * (dely1 + dely2);
732   xM[2] = xO[2] + alpha * 0.5 * (delz1 + delz2);
733 }
734 
735 /* ---------------------------------------------------------------------- */
736 
extract(const char * str,int & dim)737 void *PairLJCutTIP4PCut::extract(const char *str, int &dim)
738 {
739   dim = 0;
740   if (strcmp(str,"cut_coul") == 0) return (void *) &cut_coul;
741   dim = 2;
742   if (strcmp(str,"epsilon") == 0) return (void *) epsilon;
743   if (strcmp(str,"sigma") == 0) return (void *) sigma;
744   return nullptr;
745 }
746 /* ----------------------------------------------------------------------
747    memory usage of hneigh
748 ------------------------------------------------------------------------- */
749 
memory_usage()750 double PairLJCutTIP4PCut::memory_usage()
751 {
752   double bytes = (double)maxeatom * sizeof(double);
753   bytes += (double)maxvatom*6 * sizeof(double);
754   bytes += (double)2 * nmax * sizeof(double);
755   return bytes;
756 }
757