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 ------------------------------------------------------------------------- */
18 
19 #include "pair_lj_sf_dipole_sf_gpu.h"
20 
21 #include "atom.h"
22 #include "domain.h"
23 #include "error.h"
24 #include "force.h"
25 #include "gpu_extra.h"
26 #include "neigh_list.h"
27 #include "neigh_request.h"
28 #include "neighbor.h"
29 #include "suffix.h"
30 #include "update.h"
31 
32 #include <cmath>
33 #include <cstring>
34 
35 using namespace LAMMPS_NS;
36 
37 // External functions from cuda library for atom decomposition
38 
39 int dplsf_gpu_init(const int ntypes, double **cutsq, double **host_lj1,
40                    double **host_lj2, double **host_lj3, double **host_lj4,
41                    double *special_lj, const int nlocal,
42                    const int nall, const int max_nbors, const int maxspecial,
43                    const double cell_size, int &gpu_mode, FILE *screen,
44                    double **host_cut_ljsq, double **host_cut_coulsq,
45                    double *host_special_coul, const double qqrd2e);
46 void dplsf_gpu_clear();
47 int ** dplsf_gpu_compute_n(const int ago, const int inum, const int nall,
48                            double **host_x, int *host_type, double *sublo,
49                            double *subhi, tagint *tag, int **nspecial,
50                            tagint **special, const bool eflag, const bool vflag,
51                            const bool eatom, const bool vatom, int &host_start,
52                            int **ilist, int **jnum, const double cpu_time,
53                            bool &success, double *host_q, double **host_mu,
54                            double *boxlo, double *prd);
55 void dplsf_gpu_compute(const int ago, const int inum, const int nall,
56                        double **host_x, int *host_type, int *ilist, int *numj,
57                        int **firstneigh, const bool eflag, const bool vflag,
58                        const bool eatom, const bool vatom, int &host_start,
59                        const double cpu_time, bool &success, double *host_q,
60                        double **host_mu, const int nlocal, double *boxlo,
61                        double *prd);
62 double dplsf_gpu_bytes();
63 
64 /* ---------------------------------------------------------------------- */
65 
PairLJSFDipoleSFGPU(LAMMPS * lmp)66 PairLJSFDipoleSFGPU::PairLJSFDipoleSFGPU(LAMMPS *lmp) : PairLJSFDipoleSF(lmp),
67   gpu_mode(GPU_FORCE)
68 {
69   respa_enable = 0;
70   reinitflag = 0;
71   cpu_time = 0.0;
72   suffix_flag |= Suffix::GPU;
73   GPU_EXTRA::gpu_ready(lmp->modify, lmp->error);
74 }
75 
76 /* ----------------------------------------------------------------------
77    free all arrays
78 ------------------------------------------------------------------------- */
79 
~PairLJSFDipoleSFGPU()80 PairLJSFDipoleSFGPU::~PairLJSFDipoleSFGPU()
81 {
82   dplsf_gpu_clear();
83 }
84 
85 /* ---------------------------------------------------------------------- */
86 
compute(int eflag,int vflag)87 void PairLJSFDipoleSFGPU::compute(int eflag, int vflag)
88 {
89   ev_init(eflag,vflag);
90 
91   int nall = atom->nlocal + atom->nghost;
92   int inum, host_start;
93 
94   bool success = true;
95   int *ilist, *numneigh, **firstneigh;
96   if (gpu_mode != GPU_FORCE) {
97     double sublo[3],subhi[3];
98     if (domain->triclinic == 0) {
99       sublo[0] = domain->sublo[0];
100       sublo[1] = domain->sublo[1];
101       sublo[2] = domain->sublo[2];
102       subhi[0] = domain->subhi[0];
103       subhi[1] = domain->subhi[1];
104       subhi[2] = domain->subhi[2];
105     } else {
106       domain->bbox(domain->sublo_lamda,domain->subhi_lamda,sublo,subhi);
107     }
108     inum = atom->nlocal;
109     firstneigh = dplsf_gpu_compute_n(neighbor->ago, inum, nall, atom->x,
110                                      atom->type, sublo, subhi,
111                                      atom->tag, atom->nspecial, atom->special,
112                                      eflag, vflag, eflag_atom, vflag_atom,
113                                      host_start, &ilist, &numneigh, cpu_time,
114                                      success, atom->q, atom->mu, domain->boxlo,
115                                      domain->prd);
116   } else {
117     inum = list->inum;
118     ilist = list->ilist;
119     numneigh = list->numneigh;
120     firstneigh = list->firstneigh;
121     dplsf_gpu_compute(neighbor->ago, inum, nall, atom->x, atom->type,
122                       ilist, numneigh, firstneigh, eflag, vflag, eflag_atom,
123                       vflag_atom, host_start, cpu_time, success, atom->q,
124                       atom->mu, atom->nlocal, domain->boxlo, domain->prd);
125   }
126   if (!success)
127     error->one(FLERR,"Insufficient memory on accelerator");
128 
129   if (host_start<inum) {
130     cpu_time = MPI_Wtime();
131     cpu_compute(host_start, inum, eflag, vflag, ilist, numneigh, firstneigh);
132     cpu_time = MPI_Wtime() - cpu_time;
133   }
134 }
135 
136 /* ----------------------------------------------------------------------
137    init specific to this pair style
138 ------------------------------------------------------------------------- */
139 
init_style()140 void PairLJSFDipoleSFGPU::init_style()
141 {
142   if (!atom->q_flag || !atom->mu_flag || !atom->torque_flag)
143     error->all(FLERR,"Pair dipole/sf/gpu requires atom attributes q, mu, torque");
144 
145   if (force->newton_pair)
146     error->all(FLERR,"Pair style dipole/sf/gpu requires newton pair off");
147 
148   if (strcmp(update->unit_style,"electron") == 0)
149     error->all(FLERR,"Cannot (yet) use 'electron' units with dipoles");
150 
151   // Repeat cutsq calculation because done after call to init_style
152   double maxcut = -1.0;
153   double cut;
154   for (int i = 1; i <= atom->ntypes; i++) {
155     for (int j = i; j <= atom->ntypes; j++) {
156       if (setflag[i][j] != 0 || (setflag[i][i] != 0 && setflag[j][j] != 0)) {
157         cut = init_one(i,j);
158         cut *= cut;
159         if (cut > maxcut)
160           maxcut = cut;
161         cutsq[i][j] = cutsq[j][i] = cut;
162       } else
163         cutsq[i][j] = cutsq[j][i] = 0.0;
164     }
165   }
166   double cell_size = sqrt(maxcut) + neighbor->skin;
167 
168   int maxspecial=0;
169   if (atom->molecular != Atom::ATOMIC)
170     maxspecial=atom->maxspecial;
171   int mnf = 5e-2 * neighbor->oneatom;
172   int success = dplsf_gpu_init(atom->ntypes+1, cutsq, lj1, lj2, lj3, lj4,
173                                force->special_lj, atom->nlocal,
174                                atom->nlocal+atom->nghost, mnf, maxspecial,
175                                cell_size, gpu_mode, screen, cut_ljsq, cut_coulsq,
176                                force->special_coul, force->qqrd2e);
177   GPU_EXTRA::check_flag(success,error,world);
178 
179   if (gpu_mode == GPU_FORCE) {
180     int irequest = neighbor->request(this,instance_me);
181     neighbor->requests[irequest]->half = 0;
182     neighbor->requests[irequest]->full = 1;
183   }
184 }
185 
186 /* ---------------------------------------------------------------------- */
187 
memory_usage()188 double PairLJSFDipoleSFGPU::memory_usage()
189 {
190   double bytes = Pair::memory_usage();
191   return bytes + dplsf_gpu_bytes();
192 }
193 
194 /* ---------------------------------------------------------------------- */
195 
cpu_compute(int start,int inum,int eflag,int vflag,int * ilist,int * numneigh,int ** firstneigh)196 void PairLJSFDipoleSFGPU::cpu_compute(int start, int inum, int eflag, int vflag,
197                                   int *ilist, int *numneigh,
198                                   int **firstneigh)
199 {
200   int i,j,ii,jj,jnum,itype,jtype;
201   double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,evdwl,ecoul,fx,fy,fz;
202   double rsq,rinv,r2inv,r6inv,r3inv,r5inv;
203   double forcecoulx,forcecouly,forcecoulz,crossx,crossy,crossz;
204   double tixcoul,tiycoul,tizcoul,tjxcoul,tjycoul,tjzcoul;
205   double fq,pdotp,pidotr,pjdotr,pre1,pre2,pre3,pre4;
206   double forcelj,factor_coul,factor_lj;
207   double presf,afac,bfac,pqfac,qpfac,forceljcut,forceljsf;
208   double aforcecoulx,aforcecouly,aforcecoulz;
209   double bforcecoulx,bforcecouly,bforcecoulz;
210   double rcutlj2inv, rcutcoul2inv,rcutlj6inv;
211   int *jlist;
212 
213   evdwl = ecoul = 0.0;
214   ev_init(eflag,vflag);
215 
216   double **x = atom->x;
217   double **f = atom->f;
218   double *q = atom->q;
219   double **mu = atom->mu;
220   double **torque = atom->torque;
221   int *type = atom->type;
222   double *special_coul = force->special_coul;
223   double *special_lj = force->special_lj;
224   double qqrd2e = force->qqrd2e;
225 
226 
227   // loop over neighbors of my atoms
228 
229   for (ii = start; ii < inum; ii++) {
230     i = ilist[ii];
231     qtmp = q[i];
232     xtmp = x[i][0];
233     ytmp = x[i][1];
234     ztmp = x[i][2];
235     itype = type[i];
236     jlist = firstneigh[i];
237     jnum = numneigh[i];
238 
239     for (jj = 0; jj < jnum; jj++) {
240       j = jlist[jj];
241       factor_lj = special_lj[sbmask(j)];
242       factor_coul = special_coul[sbmask(j)];
243       j &= NEIGHMASK;
244 
245       delx = xtmp - x[j][0];
246       dely = ytmp - x[j][1];
247       delz = ztmp - x[j][2];
248       rsq = delx*delx + dely*dely + delz*delz;
249       jtype = type[j];
250 
251       if (rsq < cutsq[itype][jtype]) {
252         r2inv = 1.0/rsq;
253         rinv = sqrt(r2inv);
254 
255         // atom can have both a charge and dipole
256         // i,j = charge-charge, dipole-dipole, dipole-charge, or charge-dipole
257 
258         forcecoulx = forcecouly = forcecoulz = 0.0;
259         tixcoul = tiycoul = tizcoul = 0.0;
260         tjxcoul = tjycoul = tjzcoul = 0.0;
261 
262         if (rsq < cut_coulsq[itype][jtype]) {
263 
264           if (qtmp != 0.0 && q[j] != 0.0) {
265             pre1 = qtmp*q[j]*rinv*(r2inv-1.0/cut_coulsq[itype][jtype]);
266 
267             forcecoulx += pre1*delx;
268             forcecouly += pre1*dely;
269             forcecoulz += pre1*delz;
270           }
271 
272           if (mu[i][3] > 0.0 && mu[j][3] > 0.0) {
273             r3inv = r2inv*rinv;
274             r5inv = r3inv*r2inv;
275             rcutcoul2inv=1.0/cut_coulsq[itype][jtype];
276 
277             pdotp = mu[i][0]*mu[j][0] + mu[i][1]*mu[j][1] + mu[i][2]*mu[j][2];
278             pidotr = mu[i][0]*delx + mu[i][1]*dely + mu[i][2]*delz;
279             pjdotr = mu[j][0]*delx + mu[j][1]*dely + mu[j][2]*delz;
280 
281             afac = 1.0 - rsq*rsq * rcutcoul2inv*rcutcoul2inv;
282             pre1 = afac * ( pdotp - 3.0 * r2inv * pidotr * pjdotr );
283             aforcecoulx = pre1*delx;
284             aforcecouly = pre1*dely;
285             aforcecoulz = pre1*delz;
286 
287             bfac = 1.0 - 4.0*rsq*sqrt(rsq)*rcutcoul2inv*sqrt(rcutcoul2inv) +
288               3.0*rsq*rsq*rcutcoul2inv*rcutcoul2inv;
289             presf = 2.0 * r2inv * pidotr * pjdotr;
290             bforcecoulx = bfac * (pjdotr*mu[i][0]+pidotr*mu[j][0]-presf*delx);
291             bforcecouly = bfac * (pjdotr*mu[i][1]+pidotr*mu[j][1]-presf*dely);
292             bforcecoulz = bfac * (pjdotr*mu[i][2]+pidotr*mu[j][2]-presf*delz);
293 
294             forcecoulx += 3.0 * r5inv * ( aforcecoulx + bforcecoulx );
295             forcecouly += 3.0 * r5inv * ( aforcecouly + bforcecouly );
296             forcecoulz += 3.0 * r5inv * ( aforcecoulz + bforcecoulz );
297 
298             pre2 = 3.0 * bfac * r5inv * pjdotr;
299             pre3 = 3.0 * bfac * r5inv * pidotr;
300             pre4 = -bfac * r3inv;
301 
302             crossx = pre4 * (mu[i][1]*mu[j][2] - mu[i][2]*mu[j][1]);
303             crossy = pre4 * (mu[i][2]*mu[j][0] - mu[i][0]*mu[j][2]);
304             crossz = pre4 * (mu[i][0]*mu[j][1] - mu[i][1]*mu[j][0]);
305 
306             tixcoul += crossx + pre2 * (mu[i][1]*delz - mu[i][2]*dely);
307             tiycoul += crossy + pre2 * (mu[i][2]*delx - mu[i][0]*delz);
308             tizcoul += crossz + pre2 * (mu[i][0]*dely - mu[i][1]*delx);
309             tjxcoul += -crossx + pre3 * (mu[j][1]*delz - mu[j][2]*dely);
310             tjycoul += -crossy + pre3 * (mu[j][2]*delx - mu[j][0]*delz);
311             tjzcoul += -crossz + pre3 * (mu[j][0]*dely - mu[j][1]*delx);
312           }
313 
314           if (mu[i][3] > 0.0 && q[j] != 0.0) {
315             r3inv = r2inv*rinv;
316             r5inv = r3inv*r2inv;
317             pidotr = mu[i][0]*delx + mu[i][1]*dely + mu[i][2]*delz;
318             rcutcoul2inv=1.0/cut_coulsq[itype][jtype];
319             pre1 = 3.0 * q[j] * r5inv * pidotr * (1-rsq*rcutcoul2inv);
320             pqfac = 1.0 - 3.0*rsq*rcutcoul2inv +
321               2.0*rsq*sqrt(rsq)*rcutcoul2inv*sqrt(rcutcoul2inv);
322             pre2 = q[j] * r3inv * pqfac;
323 
324             forcecoulx += pre2*mu[i][0] - pre1*delx;
325             forcecouly += pre2*mu[i][1] - pre1*dely;
326             forcecoulz += pre2*mu[i][2] - pre1*delz;
327             tixcoul += pre2 * (mu[i][1]*delz - mu[i][2]*dely);
328             tiycoul += pre2 * (mu[i][2]*delx - mu[i][0]*delz);
329             tizcoul += pre2 * (mu[i][0]*dely - mu[i][1]*delx);
330           }
331 
332           if (mu[j][3] > 0.0 && qtmp != 0.0) {
333             r3inv = r2inv*rinv;
334             r5inv = r3inv*r2inv;
335             pjdotr = mu[j][0]*delx + mu[j][1]*dely + mu[j][2]*delz;
336             rcutcoul2inv=1.0/cut_coulsq[itype][jtype];
337             pre1 = 3.0 * qtmp * r5inv * pjdotr * (1-rsq*rcutcoul2inv);
338             qpfac = 1.0 - 3.0*rsq*rcutcoul2inv +
339               2.0*rsq*sqrt(rsq)*rcutcoul2inv*sqrt(rcutcoul2inv);
340             pre2 = qtmp * r3inv * qpfac;
341 
342             forcecoulx += pre1*delx - pre2*mu[j][0];
343             forcecouly += pre1*dely - pre2*mu[j][1];
344             forcecoulz += pre1*delz - pre2*mu[j][2];
345             tjxcoul += -pre2 * (mu[j][1]*delz - mu[j][2]*dely);
346             tjycoul += -pre2 * (mu[j][2]*delx - mu[j][0]*delz);
347             tjzcoul += -pre2 * (mu[j][0]*dely - mu[j][1]*delx);
348           }
349         }
350 
351         // LJ interaction
352 
353         if (rsq < cut_ljsq[itype][jtype]) {
354           r6inv = r2inv*r2inv*r2inv;
355           forceljcut = r6inv*(lj1[itype][jtype]*r6inv-lj2[itype][jtype])*r2inv;
356 
357           rcutlj2inv = 1.0 / cut_ljsq[itype][jtype];
358           rcutlj6inv = rcutlj2inv * rcutlj2inv * rcutlj2inv;
359           forceljsf = (lj1[itype][jtype]*rcutlj6inv - lj2[itype][jtype]) *
360           rcutlj6inv * rcutlj2inv;
361 
362           forcelj = factor_lj * (forceljcut - forceljsf);
363         } else forcelj = 0.0;
364 
365         // total force
366 
367         fq = factor_coul*qqrd2e;
368         fx = fq*forcecoulx + delx*forcelj;
369         fy = fq*forcecouly + dely*forcelj;
370         fz = fq*forcecoulz + delz*forcelj;
371 
372         // force & torque accumulation
373 
374         f[i][0] += fx;
375         f[i][1] += fy;
376         f[i][2] += fz;
377         torque[i][0] += fq*tixcoul;
378         torque[i][1] += fq*tiycoul;
379         torque[i][2] += fq*tizcoul;
380 
381         if (eflag) {
382           if (rsq < cut_coulsq[itype][jtype]) {
383             ecoul = qtmp*q[j]*rinv*
384               pow((1.0-sqrt(rsq)/sqrt(cut_coulsq[itype][jtype])),2);
385             if (mu[i][3] > 0.0 && mu[j][3] > 0.0)
386               ecoul += bfac * (r3inv*pdotp - 3.0*r5inv*pidotr*pjdotr);
387             if (mu[i][3] > 0.0 && q[j] != 0.0)
388               ecoul += -q[j]*r3inv * pqfac * pidotr;
389             if (mu[j][3] > 0.0 && qtmp != 0.0)
390               ecoul += qtmp*r3inv * qpfac * pjdotr;
391             ecoul *= factor_coul*qqrd2e;
392           } else ecoul = 0.0;
393 
394           if (rsq < cut_ljsq[itype][jtype]) {
395             evdwl = r6inv*(lj3[itype][jtype]*r6inv-lj4[itype][jtype]) +
396               rcutlj6inv*(6*lj3[itype][jtype]*rcutlj6inv-3*lj4[itype][jtype])*
397               rsq*rcutlj2inv +
398               rcutlj6inv*(-7*lj3[itype][jtype]*rcutlj6inv+4*lj4[itype][jtype]);
399             evdwl *= factor_lj;
400           } else evdwl = 0.0;
401         }
402 
403         if (evflag) ev_tally_xyz_full(i,evdwl,ecoul,
404                                       fx,fy,fz,delx,dely,delz);
405       }
406     }
407   }
408 }
409