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: Zhen Li (Brown University)
17 Email: zhen_li@brown.edu
18 ------------------------------------------------------------------------- */
19
20 #include "pair_mdpd.h"
21
22 #include <cmath>
23 #include <ctime>
24 #include "atom.h"
25 #include "comm.h"
26 #include "update.h"
27 #include "force.h"
28 #include "neighbor.h"
29 #include "neigh_list.h"
30 #include "random_mars.h"
31 #include "citeme.h"
32 #include "memory.h"
33 #include "error.h"
34
35
36 using namespace LAMMPS_NS;
37
38 #define EPSILON 1.0e-10
39
40 static const char cite_pair_mdpd[] =
41 "pair mdpd command:\n\n"
42 "@Article{ZLi2013_POF,\n"
43 " author = {Li, Z. and Hu, G.H. and Wang, Z.L. and Ma Y.B. and Zhou, Z.W.},\n"
44 " title = {Three dimensional flow structures in a moving droplet on substrate: a dissipative particle dynamics study},\n"
45 " journal = {Physics of Fluids},\n"
46 " year = {2013},\n"
47 " volume = {25},\n"
48 " pages = {072103}\n"
49 "}\n\n";
50
51 /* ---------------------------------------------------------------------- */
52
PairMDPD(LAMMPS * lmp)53 PairMDPD::PairMDPD(LAMMPS *lmp) : Pair(lmp)
54 {
55 if (lmp->citeme) lmp->citeme->add(cite_pair_mdpd);
56
57 writedata = 1;
58 random = nullptr;
59 }
60
61 /* ---------------------------------------------------------------------- */
62
~PairMDPD()63 PairMDPD::~PairMDPD()
64 {
65 if (allocated) {
66 memory->destroy(setflag);
67 memory->destroy(cutsq);
68
69 memory->destroy(cut);
70 memory->destroy(cut_r);
71 memory->destroy(A_att);
72 memory->destroy(B_rep);
73 memory->destroy(gamma);
74 memory->destroy(sigma);
75 }
76 if (random) delete random;
77 }
78
79 /* ---------------------------------------------------------------------- */
80
compute(int eflag,int vflag)81 void PairMDPD::compute(int eflag, int vflag)
82 {
83 int i,j,ii,jj,inum,jnum,itype,jtype;
84 double xtmp,ytmp,ztmp,delx,dely,delz,evdwl,fpair;
85 double vxtmp,vytmp,vztmp,delvx,delvy,delvz;
86 double rsq,r,rinv,dot,wc,wc_r, wr,randnum,factor_dpd;
87 int *ilist,*jlist,*numneigh,**firstneigh;
88 double rhoi, rhoj;
89
90 evdwl = 0.0;
91 ev_init(eflag,vflag);
92
93 double **x = atom->x;
94 double **v = atom->v;
95 double **f = atom->f;
96 double *rho= atom->rho;
97 int *type = atom->type;
98 int nlocal = atom->nlocal;
99 double *special_lj = force->special_lj;
100 int newton_pair = force->newton_pair;
101 double dtinvsqrt = 1.0/sqrt(update->dt);
102
103 inum = list->inum;
104 ilist = list->ilist;
105 numneigh = list->numneigh;
106 firstneigh = list->firstneigh;
107
108 // loop over neighbors of my atoms
109
110 for (ii = 0; ii < inum; ii++) {
111 i = ilist[ii];
112 xtmp = x[i][0];
113 ytmp = x[i][1];
114 ztmp = x[i][2];
115 vxtmp = v[i][0];
116 vytmp = v[i][1];
117 vztmp = v[i][2];
118 itype = type[i];
119 jlist = firstneigh[i];
120 jnum = numneigh[i];
121 rhoi = rho[i];
122 for (jj = 0; jj < jnum; jj++) {
123 j = jlist[jj];
124 factor_dpd = special_lj[sbmask(j)];
125 j &= NEIGHMASK;
126
127 delx = xtmp - x[j][0];
128 dely = ytmp - x[j][1];
129 delz = ztmp - x[j][2];
130 rsq = delx*delx + dely*dely + delz*delz;
131 jtype = type[j];
132
133 if (rsq < cutsq[itype][jtype]) {
134 r = sqrt(rsq);
135 if (r < EPSILON) continue; // r can be 0.0 in MDPD systems
136 rinv = 1.0/r;
137 delvx = vxtmp - v[j][0];
138 delvy = vytmp - v[j][1];
139 delvz = vztmp - v[j][2];
140 dot = delx*delvx + dely*delvy + delz*delvz;
141
142 wc = 1.0 - r/cut[itype][jtype];
143 wc_r = 1.0 - r/cut_r[itype][jtype];
144 wc_r = MAX(wc_r,0.0);
145 wr = wc;
146
147 rhoj = rho[j];
148 randnum = random->gaussian();
149
150 // conservative force = A_att * wc + B_rep*(rhoi+rhoj)*wc_r
151 // drag force = -gamma * wr^2 * (delx dot delv) / r
152 // random force = sigma * wr * rnd * dtinvsqrt;
153
154 fpair = A_att[itype][jtype]*wc + B_rep[itype][jtype]*(rhoi+rhoj)*wc_r;
155 fpair -= gamma[itype][jtype]*wr*wr*dot*rinv;
156 fpair += sigma[itype][jtype]*wr*randnum*dtinvsqrt;
157 fpair *= factor_dpd*rinv;
158
159 f[i][0] += delx*fpair;
160 f[i][1] += dely*fpair;
161 f[i][2] += delz*fpair;
162 if (newton_pair || j < nlocal) {
163 f[j][0] -= delx*fpair;
164 f[j][1] -= dely*fpair;
165 f[j][2] -= delz*fpair;
166 }
167
168 if (eflag) {
169 // unshifted eng of conservative term:
170 // eng shifted to 0.0 at cutoff
171 evdwl = 0.5*A_att[itype][jtype]*cut[itype][jtype] * wr*wr + 0.5*B_rep[itype][jtype]*cut_r[itype][jtype]*(rhoi+rhoj)*wc_r*wc_r;
172 evdwl *= factor_dpd;
173 }
174
175 if (evflag) ev_tally(i,j,nlocal,newton_pair,evdwl,0.0,fpair,delx,dely,delz);
176 }
177 }
178 }
179 if (vflag_fdotr) virial_fdotr_compute();
180 }
181
182 /* ----------------------------------------------------------------------
183 allocate all arrays
184 ------------------------------------------------------------------------- */
185
allocate()186 void PairMDPD::allocate()
187 {
188 int i,j;
189 allocated = 1;
190 int n = atom->ntypes;
191
192 memory->create(setflag,n+1,n+1,"pair:setflag");
193 for (i = 1; i <= n; i++)
194 for (j = i; j <= n; j++)
195 setflag[i][j] = 0;
196
197 memory->create(cutsq,n+1,n+1,"pair:cutsq");
198
199 memory->create(cut,n+1,n+1,"pair:cut");
200 memory->create(cut_r,n+1,n+1,"pair:cut_r");
201 memory->create(A_att,n+1,n+1,"pair:A_att");
202 memory->create(B_rep,n+1,n+1,"pair:B_rep");
203 memory->create(gamma,n+1,n+1,"pair:gamma");
204 memory->create(sigma,n+1,n+1,"pair:sigma");
205 }
206
207 /* ----------------------------------------------------------------------
208 global settings
209 ------------------------------------------------------------------------- */
210
settings(int narg,char ** arg)211 void PairMDPD::settings(int narg, char **arg)
212 {
213 if (narg != 3) error->all(FLERR,"Illegal pair_style command");
214
215 temperature = utils::numeric(FLERR,arg[0],false,lmp);
216 cut_global = utils::numeric(FLERR,arg[1],false,lmp);
217 seed = utils::inumeric(FLERR,arg[2],false,lmp);
218
219 // initialize Marsaglia RNG with processor-unique seed
220
221 if (seed <= 0) {
222 struct timespec time;
223 clock_gettime( CLOCK_REALTIME, &time );
224 seed = time.tv_nsec; // if seed is non-positive, get the current time as the seed
225 }
226 delete random;
227 random = new RanMars(lmp,(seed + comm->me) % 900000000);
228
229 // reset cutoffs that have been explicitly set
230
231 if (allocated) {
232 int i,j;
233 for (i = 1; i <= atom->ntypes; i++)
234 for (j = i+1; j <= atom->ntypes; j++)
235 if (setflag[i][j]) cut[i][j] = cut_global;
236 }
237 }
238
239 /* ----------------------------------------------------------------------
240 set coeffs for one or more type pairs
241 ------------------------------------------------------------------------- */
242
coeff(int narg,char ** arg)243 void PairMDPD::coeff(int narg, char **arg)
244 {
245 if (narg != 7 ) error->all(FLERR,"Incorrect args for pair coefficients\n itype jtype A B gamma cutA cutB");
246 if (!allocated) allocate();
247
248 int ilo,ihi,jlo,jhi;
249 utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error);
250 utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error);
251
252 double A_one = utils::numeric(FLERR,arg[2],false,lmp);
253 double B_one = utils::numeric(FLERR,arg[3],false,lmp);
254 double gamma_one = utils::numeric(FLERR,arg[4],false,lmp);
255 double cut_one = utils::numeric(FLERR,arg[5],false,lmp);
256 double cut_two = utils::numeric(FLERR,arg[6],false,lmp);
257
258 if (cut_one < cut_two) error->all(FLERR,"Incorrect args for pair coefficients\n cutA should be larger than cutB.");
259
260 int count = 0;
261 for (int i = ilo; i <= ihi; i++) {
262 for (int j = MAX(jlo,i); j <= jhi; j++) {
263 A_att[i][j] = A_one;
264 B_rep[i][j] = B_one;
265 gamma[i][j] = gamma_one;
266 cut[i][j] = cut_one;
267 cut_r[i][j] = cut_two;
268 setflag[i][j] = 1;
269 count++;
270 }
271 }
272 if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients");
273 }
274
275 /* ----------------------------------------------------------------------
276 init specific to this pair style
277 ------------------------------------------------------------------------- */
278
init_style()279 void PairMDPD::init_style()
280 {
281 if (comm->ghost_velocity == 0)
282 error->all(FLERR,"Pair mdpd requires ghost atoms store velocity");
283
284 if (!atom->rho_flag)
285 error->all(FLERR,"Pair style mdpd requires atom attribute rho");
286
287 // if newton off, forces between atoms ij will be double computed
288 // using different random numbers
289
290 if (force->newton_pair == 0 && comm->me == 0) error->warning(FLERR,
291 "Pair mdpd needs newton pair on for momentum conservation");
292
293 neighbor->request(this,instance_me);
294 }
295
296 /* ----------------------------------------------------------------------
297 init for one type pair i,j and corresponding j,i
298 ------------------------------------------------------------------------- */
299
init_one(int i,int j)300 double PairMDPD::init_one(int i, int j)
301 {
302 if (setflag[i][j] == 0) error->all(FLERR,"All pair coeffs are not set");
303
304 sigma[i][j] = sqrt(2.0*force->boltz*temperature*gamma[i][j]);
305
306 cut[j][i] = cut[i][j];
307 cut_r[j][i] = cut_r[i][j];
308 A_att[j][i] = A_att[i][j];
309 B_rep[j][i] = B_rep[i][j];
310 gamma[j][i] = gamma[i][j];
311 sigma[j][i] = sigma[i][j];
312
313 return cut[i][j];
314 }
315
316 /* ----------------------------------------------------------------------
317 proc 0 writes to restart file
318 ------------------------------------------------------------------------- */
319
write_restart(FILE * fp)320 void PairMDPD::write_restart(FILE *fp)
321 {
322 write_restart_settings(fp);
323
324 int i,j;
325 for (i = 1; i <= atom->ntypes; i++)
326 for (j = i; j <= atom->ntypes; j++) {
327 fwrite(&setflag[i][j],sizeof(int),1,fp);
328 if (setflag[i][j]) {
329 fwrite(&A_att[i][j],sizeof(double),1,fp);
330 fwrite(&B_rep[i][j],sizeof(double),1,fp);
331 fwrite(&gamma[i][j],sizeof(double),1,fp);
332 fwrite(&cut[i][j],sizeof(double),1,fp);
333 fwrite(&cut_r[i][j],sizeof(double),1,fp);
334 }
335 }
336 }
337
338 /* ----------------------------------------------------------------------
339 proc 0 reads from restart file, bcasts
340 ------------------------------------------------------------------------- */
341
read_restart(FILE * fp)342 void PairMDPD::read_restart(FILE *fp)
343 {
344 read_restart_settings(fp);
345
346 allocate();
347
348 int i,j;
349 int me = comm->me;
350 for (i = 1; i <= atom->ntypes; i++)
351 for (j = i; j <= atom->ntypes; j++) {
352 if (me == 0) utils::sfread(FLERR,&setflag[i][j],sizeof(int),1,fp,nullptr,error);
353 MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world);
354 if (setflag[i][j]) {
355 if (me == 0) {
356 utils::sfread(FLERR,&A_att[i][j],sizeof(double),1,fp,nullptr,error);
357 utils::sfread(FLERR,&B_rep[i][j],sizeof(double),1,fp,nullptr,error);
358 utils::sfread(FLERR,&gamma[i][j],sizeof(double),1,fp,nullptr,error);
359 utils::sfread(FLERR,&cut[i][j],sizeof(double),1,fp,nullptr,error);
360 utils::sfread(FLERR,&cut_r[i][j],sizeof(double),1,fp,nullptr,error);
361 }
362 MPI_Bcast(&A_att[i][j],1,MPI_DOUBLE,0,world);
363 MPI_Bcast(&B_rep[i][j],1,MPI_DOUBLE,0,world);
364 MPI_Bcast(&gamma[i][j],1,MPI_DOUBLE,0,world);
365 MPI_Bcast(&cut[i][j],1,MPI_DOUBLE,0,world);
366 MPI_Bcast(&cut_r[i][j],1,MPI_DOUBLE,0,world);
367 }
368 }
369 }
370
371 /* ----------------------------------------------------------------------
372 proc 0 writes to restart file
373 ------------------------------------------------------------------------- */
374
write_restart_settings(FILE * fp)375 void PairMDPD::write_restart_settings(FILE *fp)
376 {
377 fwrite(&temperature,sizeof(double),1,fp);
378 fwrite(&cut_global,sizeof(double),1,fp);
379 fwrite(&seed,sizeof(int),1,fp);
380 fwrite(&mix_flag,sizeof(int),1,fp);
381 }
382
383 /* ----------------------------------------------------------------------
384 proc 0 reads from restart file, bcasts
385 ------------------------------------------------------------------------- */
386
read_restart_settings(FILE * fp)387 void PairMDPD::read_restart_settings(FILE *fp)
388 {
389 if (comm->me == 0) {
390 utils::sfread(FLERR,&temperature,sizeof(double),1,fp,nullptr,error);
391 utils::sfread(FLERR,&cut_global,sizeof(double),1,fp,nullptr,error);
392 utils::sfread(FLERR,&seed,sizeof(int),1,fp,nullptr,error);
393 utils::sfread(FLERR,&mix_flag,sizeof(int),1,fp,nullptr,error);
394 }
395 MPI_Bcast(&temperature,1,MPI_DOUBLE,0,world);
396 MPI_Bcast(&cut_global,1,MPI_DOUBLE,0,world);
397 MPI_Bcast(&seed,1,MPI_INT,0,world);
398 MPI_Bcast(&mix_flag,1,MPI_INT,0,world);
399
400 // initialize Marsaglia RNG with processor-unique seed
401 // same seed that pair_style command initially specified
402
403 if (random) delete random;
404 random = new RanMars(lmp,seed + comm->me);
405 }
406
407 /* ----------------------------------------------------------------------
408 proc 0 writes to data file
409 ------------------------------------------------------------------------- */
410
write_data(FILE * fp)411 void PairMDPD::write_data(FILE *fp)
412 {
413 for (int i = 1; i <= atom->ntypes; i++)
414 fprintf(fp,"%d %g %g %g\n",i,A_att[i][i],B_rep[i][i],gamma[i][i]);
415 }
416
417 /* ----------------------------------------------------------------------
418 proc 0 writes all pairs to data file
419 ------------------------------------------------------------------------- */
420
write_data_all(FILE * fp)421 void PairMDPD::write_data_all(FILE *fp)
422 {
423 for (int i = 1; i <= atom->ntypes; i++)
424 for (int j = i; j <= atom->ntypes; j++)
425 fprintf(fp,"%d %d %g %g %g %g %g\n",i,j,A_att[i][j],B_rep[i][j],gamma[i][j],cut[i][j],cut_r[i][j]);
426 }
427
428