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 #include "compute_pressure_cylinder.h"
16
17 #include <cmath>
18 #include "atom.h"
19 #include "update.h"
20 #include "force.h"
21 #include "pair.h"
22 #include "neighbor.h"
23 #include "neigh_request.h"
24 #include "neigh_list.h"
25 #include "memory.h"
26 #include "error.h"
27 #include "citeme.h"
28 #include "domain.h"
29 #include "math_const.h"
30
31 using namespace LAMMPS_NS;
32 using namespace MathConst;
33
34 static const char cite_compute_pressure_cylinder[] =
35 "compute pressure/cylinder:\n\n"
36 "@Article{Addington,\n"
37 " author = {C. K. Addington, Y. Long, K. E. Gubbins},\n"
38 " title = {The pressure in interfaces having cylindrical geometry},\n"
39 " journal = {J.~Chem.~Phys.},\n"
40 " year = 2018,\n"
41 " volume = 149,\n"
42 " pages = {084109}\n"
43 "}\n\n";
44
45 /*+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
46 Calculate the configurational components of the pressure tensor in
47 cylindrical geometry, according to the formulation of Addington et al. (2018)
48 +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
49
ComputePressureCyl(LAMMPS * lmp,int narg,char ** arg)50 ComputePressureCyl::ComputePressureCyl(LAMMPS *lmp, int narg, char **arg) :
51 Compute(lmp, narg, arg),
52 Pr_temp(nullptr), Pr_all(nullptr), Pz_temp(nullptr), Pz_all(nullptr), Pphi_temp(nullptr),
53 Pphi_all(nullptr), R(nullptr), Rinv(nullptr), R2(nullptr), PrAinv(nullptr), PzAinv(nullptr),
54 R2kin(nullptr), density_temp(nullptr), invVbin(nullptr), density_all(nullptr),
55 tangent(nullptr), ephi_x(nullptr), ephi_y(nullptr), binz(nullptr)
56 {
57 if (lmp->citeme) lmp->citeme->add(cite_compute_pressure_cylinder);
58 if (narg != 7) error->all(FLERR,"Illegal compute pressure/cylinder command");
59
60 zlo=utils::numeric(FLERR,arg[3],false,lmp);
61 zhi=utils::numeric(FLERR,arg[4],false,lmp);
62 Rmax=utils::numeric(FLERR,arg[5],false,lmp);
63 bin_width=utils::numeric(FLERR,arg[6],false,lmp);
64
65 if ((bin_width <= 0.0) || (bin_width > Rmax))
66 error->all(FLERR,"Illegal compute pressure/cylinder command");
67 if ((zhi < zlo) || ((zhi-zlo) < bin_width))
68 error->all(FLERR,"Illegal compute pressure/cylinder command");
69 if ((zhi > domain->boxhi[2]) || (zlo < domain->boxlo[2]))
70 error->all(FLERR,"Illegal compute pressure/cylinder command");
71
72 nbins=(int)(Rmax/bin_width);
73 nzbins=(int)((zhi-zlo)/bin_width);
74
75 // NOTE: at 2^22 = 4.2M bins, we will be close to exhausting allocatable
76 // memory on a 32-bit environment. so we use this as an upper limit.
77
78 if ((nbins < 1) || (nzbins < 1) || (nbins > 2<<22) || (nzbins > 2<<22))
79 error->all(FLERR,"Illegal compute pressure/cylinder command");
80
81 array_flag=1;
82 vector_flag=0;
83 extarray=0;
84 size_array_cols = 5; // r, number density, Pr, Pphi, Pz
85 size_array_rows = nbins;
86
87 Pr_temp = new double[nbins];
88 Pr_all = new double[nbins];
89 Pz_temp = new double[nbins];
90 Pz_all = new double[nbins];
91 Pphi_temp = new double[nbins];
92 Pphi_all = new double[nbins];
93 R = new double[nbins];
94 R2 = new double[nbins];
95 PrAinv = new double[nbins];
96 PzAinv = new double[nbins];
97 Rinv = new double[nbins];
98 binz = new double[nzbins];
99
100 R2kin = new double[nbins];
101 density_temp = new double[nbins];
102 invVbin = new double[nbins];
103 density_all = new double[nbins];
104
105 memory->create(array,nbins,5,"PN:array");
106
107 nphi=360;
108 tangent = new double[nphi];
109 ephi_x = new double[nphi];
110 ephi_y = new double[nphi];
111
112 nktv2p = force->nktv2p;
113
114 }
115
116 /* ---------------------------------------------------------------------- */
117
~ComputePressureCyl()118 ComputePressureCyl::~ComputePressureCyl()
119 {
120 // count all of these for memory usage
121 memory->destroy(array);
122 delete [] R;
123 delete [] Rinv;
124 delete [] R2;
125 delete [] R2kin;
126 delete [] invVbin;
127 delete [] density_temp;
128 delete [] density_all;
129 delete [] tangent;
130 delete [] ephi_x;
131 delete [] ephi_y;
132 delete [] Pr_temp;
133 delete [] Pr_all;
134 delete [] Pz_temp;
135 delete [] Pz_all;
136 delete [] Pphi_temp;
137 delete [] Pphi_all;
138 delete [] PrAinv;
139 delete [] PzAinv;
140 delete [] binz;
141 }
142
143 /* ---------------------------------------------------------------------- */
144
init()145 void ComputePressureCyl::init()
146 {
147 if (force->pair == nullptr)
148 error->all(FLERR,"No pair style is defined for compute pressure/cylinder");
149 if (force->pair->single_enable == 0)
150 error->all(FLERR,"Pair style does not support compute pressure/cylinder");
151
152 double phi;
153
154 for (int iphi = 0; iphi < nphi; iphi++) {
155 phi=((double)iphi)*MY_PI/180.0;
156 tangent[iphi]=tan(phi);
157 ephi_x[iphi]=-sin(phi);
158 ephi_y[iphi]=cos(phi);
159 }
160
161 for (int iq = 0; iq < nbins; iq++) {
162 R[iq]=((double)iq+0.5)*bin_width;
163 Rinv[iq]=1.0/R[iq];
164 R2[iq]=R[iq]*R[iq];
165 R2kin[iq]=(((double)iq)+1.0)*bin_width;
166 R2kin[iq]*=R2kin[iq];
167 PrAinv[iq]=1.0/(2.0*MY_PI*(zhi-zlo)*R[iq]);
168 }
169 PphiAinv=1.0/((zhi-zlo)*bin_width*2.0*(double)nphi);
170
171 invVbin[0]=1.0/((zhi-zlo)*MY_PI*R2kin[0]);
172 PzAinv[0]=1.0/(MY_PI*R2kin[0]*((double)nzbins));
173
174 for (int jq = 1; jq < nbins; jq++) {
175 invVbin[jq]=1.0/((zhi-zlo)*MY_PI*(R2kin[jq]-R2kin[jq-1]));
176 PzAinv[jq]=1.0/(MY_PI*(R2kin[jq]-R2kin[jq-1])*((double)nzbins));
177 }
178
179 // need an occasional half neighbor list
180 int irequest = neighbor->request(this,instance_me);
181 neighbor->requests[irequest]->pair = 0;
182 neighbor->requests[irequest]->compute = 1;
183 neighbor->requests[irequest]->occasional = 1;
184
185 for (int zzz = 0; zzz < nzbins; zzz++) binz[zzz]=(((double)zzz)+0.5)*bin_width+zlo;
186
187 }
188
189 /* ---------------------------------------------------------------------- */
190
init_list(int,NeighList * ptr)191 void ComputePressureCyl::init_list(int /* id */, NeighList *ptr)
192 {
193 list = ptr;
194 }
195
196 /* ---------------------------------------------------------------------- */
197
198
199 /* ----------------------------------------------------------------------
200 count pairs and compute pair info on this proc
201 only count pair once if newton_pair is off
202 both atom I,J must be in group
203 if flag is set, compute requested info about pair
204 ------------------------------------------------------------------------- */
205
compute_array()206 void ComputePressureCyl::compute_array()
207 {
208 invoked_array = update->ntimestep;
209
210 int ibin;
211
212 // clear pressures
213 for (ibin = 0; ibin < nbins; ibin++) {
214 density_temp[ibin]=0.0;
215 density_all[ibin]=0.0;
216 Pr_temp[ibin]=0.0;
217 Pr_all[ibin]=0.0;
218 Pphi_temp[ibin]=0.0;
219 Pphi_all[ibin]=0.0;
220 Pz_temp[ibin]=0.0;
221 Pz_all[ibin]=0.0;
222 }
223
224 // what processor am I?
225 int me;
226 MPI_Comm_rank(world,&me);
227
228 int i,j,ii,jj,inum,jnum,itype,jtype;
229 tagint itag,jtag;
230 double xtmp,ytmp,ztmp,delx,dely,delz;
231 double rsq,fpair,factor_coul,factor_lj;
232 int *ilist,*jlist,*numneigh,**firstneigh;
233
234 double **x = atom->x;
235 tagint *tag = atom->tag;
236 int *type = atom->type;
237 int *mask = atom->mask;
238 int nlocal = atom->nlocal;
239 double *special_coul = force->special_coul;
240 double *special_lj = force->special_lj;
241 int newton_pair = force->newton_pair;
242
243 // invoke half neighbor list (will copy or build if necessary)
244 neighbor->build_one(list);
245
246 inum = list->inum;
247 ilist = list->ilist;
248 numneigh = list->numneigh;
249 firstneigh = list->firstneigh;
250
251 // calculate number density (by radius)
252 double temp_R2;
253 for (i = 0; i < nlocal; i++) if ((x[i][2] < zhi) && (x[i][2] > zlo)) {
254 temp_R2=x[i][0]*x[i][0]+x[i][1]*x[i][1];
255 if (temp_R2 > R2kin[nbins-1]) continue; // outside of Rmax
256
257 for (j = 0; j < nbins; j++) if (temp_R2 < R2kin[j]) break;
258
259 density_temp[j]+=invVbin[j];
260 }
261 MPI_Allreduce(density_temp,density_all,nbins,MPI_DOUBLE,MPI_SUM,world);
262 for (i = 0; i < nbins; i++) array[i][1]=density_all[i]; // NEW
263
264 // loop over neighbors of my atoms
265 // skip if I or J are not in group
266 // for newton = 0 and J = ghost atom,
267 // need to insure I,J pair is only output by one proc
268 // use same itag,jtag logic as in Neighbor::neigh_half_nsq()
269 // for flag = 0, just count pair interactions within force cutoff
270 // for flag = 1, calculate requested output fields
271
272 Pair *pair = force->pair;
273 double **cutsq = force->pair->cutsq;
274
275 double r1=0.0;
276 double r2=0.0;
277 double risq,rjsq;
278 double A,B,C,D;
279 double alpha1,alpha2;
280 double xi,yi,zi,dx,dy,dz;
281 double xR,yR,zR,fn;
282 double alpha,xL,yL,zL,L2,ftphi,ftz;
283 double sqrtD;
284
285 for (ii = 0; ii < inum; ii++) {
286 i = ilist[ii];
287 if (!(mask[i] & groupbit)) continue;
288
289 xtmp = x[i][0];
290 ytmp = x[i][1];
291 ztmp = x[i][2];
292 itag = tag[i];
293 itype = type[i];
294 jlist = firstneigh[i];
295 jnum = numneigh[i];
296
297 r1=x[i][0]*x[i][0]+x[i][1]*x[i][1];
298
299 for (jj = 0; jj < jnum; jj++) {
300 j = jlist[jj];
301 factor_lj = special_lj[sbmask(j)];
302 factor_coul = special_coul[sbmask(j)];
303 j &= NEIGHMASK;
304
305 if (!(mask[j] & groupbit)) continue;
306
307 // itag = jtag is possible for long cutoffs that include images of self
308 // do calculation only on appropriate processor
309 if (newton_pair == 0 && j >= nlocal) {
310 jtag = tag[j];
311 if (itag > jtag) {
312 if ((itag+jtag) % 2 == 0) continue;
313 } else if (itag < jtag) {
314 if ((itag+jtag) % 2 == 1) continue;
315 } else {
316 if (x[j][2] < ztmp) continue;
317 if (x[j][2] == ztmp) {
318 if (x[j][1] < ytmp) continue;
319 if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
320 }
321 }
322 }
323
324 delx = xtmp - x[j][0];
325 dely = ytmp - x[j][1];
326 delz = ztmp - x[j][2];
327
328 r2=x[j][0]*x[j][0]+x[j][1]*x[j][1];
329
330 // ri is smaller of r1 and r2
331 if (r2 < r1) {
332 risq=r2;
333 rjsq=r1;
334 xi=x[j][0];
335 yi=x[j][1];
336 zi=x[j][2];
337 dx=x[i][0]-x[j][0];
338 dy=x[i][1]-x[j][1];
339 dz=x[i][2]-x[j][2];
340 } else {
341 risq=r1;
342 rjsq=r2;
343 xi=x[i][0];
344 yi=x[i][1];
345 zi=x[i][2];
346 dx=x[j][0]-x[i][0];
347 dy=x[j][1]-x[i][1];
348 dz=x[j][2]-x[i][2];
349 }
350
351 rsq = delx*delx + dely*dely + delz*delz;
352 jtype = type[j];
353 if (rsq >= cutsq[itype][jtype]) continue;
354
355 pair->single(i,j,itype,jtype,rsq,factor_coul,factor_lj,fpair);
356
357 A=dx*dx+dy*dy;
358 B=2.0*(xi*dx+yi*dy);
359
360 // normal pressure contribution P_rhorho
361 for (ibin = 0; ibin < nbins; ibin++) {
362 // completely inside of R
363 if (rjsq < R2[ibin]) continue;
364
365 C=risq-R2[ibin];
366 D=B*B-4.0*A*C;
367
368 // completely outside of R
369 if (D < 0.0) continue;
370
371 sqrtD=sqrt(D);
372 alpha1=0.5*(-B+sqrtD)/A;
373 alpha2=0.5*(-B-sqrtD)/A;
374
375 if ((alpha1 > 0.0) && (alpha1 < 1.0)) {
376 zR=zi+alpha1*dz;
377 if ((zR < zhi) && (zR > zlo))
378 {
379 xR=xi+alpha1*dx;
380 yR=yi+alpha1*dy;
381 fn=fpair*fabs(xR*dx+yR*dy);
382
383 Pr_temp[ibin]+=fn;
384 }
385 }
386 if ((alpha2 > 0.0) && (alpha2 < 1.0)) {
387 zR=zi+alpha2*dz;
388 if ((zR < zhi) && (zR > zlo)) {
389 xR=xi+alpha2*dx;
390 yR=yi+alpha2*dy;
391 fn=fpair*fabs(xR*dx+yR*dy);
392
393 Pr_temp[ibin]+=fn;
394 }
395 }
396 }
397
398 // azimuthal pressure contribution (P_phiphi)
399 for (int iphi = 0; iphi < nphi; iphi++) {
400 alpha=(yi-xi*tangent[iphi])/(dx*tangent[iphi]-dy);
401
402 // no intersection with phi surface
403 if ((alpha >= 1.0) || (alpha <= 0.0)) continue;
404
405 // no contribution (outside of averaging region)
406 zL=zi+alpha*dz;
407 if ((zL > zhi) || (zL < zlo)) continue;
408
409 xL=xi+alpha*dx;
410 yL=yi+alpha*dy;
411
412 L2=xL*xL+yL*yL;
413
414 // no intersection (outside of Rmax)
415 if (L2 > R2kin[nbins-1]) continue;
416
417 ftphi=fabs(dx*ephi_x[iphi]+dy*ephi_y[iphi])*fpair;
418
419 // add to appropriate bin
420 for (ibin = 0; ibin < nbins; ibin++) if (L2 < R2kin[ibin]) {
421 Pphi_temp[ibin]+=ftphi;
422 break;
423 }
424 }
425
426 // z pressure contribution (P_zz)
427 for (int zbin = 0; zbin < nzbins; zbin++) {
428 // check if interaction contributes
429 if ((x[i][2] > binz[zbin]) && (x[j][2] > binz[zbin])) continue;
430 if ((x[i][2] < binz[zbin]) && (x[j][2] < binz[zbin])) continue;
431
432 alpha=(binz[zbin]-zi)/dz;
433
434 xL=xi+alpha*dx;
435 yL=yi+alpha*dy;
436
437 L2=xL*xL+yL*yL;
438
439 if (L2 > R2kin[nbins-1]) continue;
440
441 ftz=fabs(dz)*fpair;
442
443 // add to appropriate bin
444 for (ibin = 0; ibin < nbins; ibin++) if (L2 < R2kin[ibin]) {
445 Pz_temp[ibin]+=ftz;
446 break;
447 }
448 }
449 }
450 }
451
452 // calculate pressure (force over area)
453 for (ibin = 0; ibin < nbins; ibin++) {
454 Pr_temp[ibin]*=PrAinv[ibin]*Rinv[ibin];
455 Pphi_temp[ibin]*=PphiAinv;
456 Pz_temp[ibin]*=PzAinv[ibin];
457 }
458
459 // communicate these values across processors
460 MPI_Allreduce(Pr_temp,Pr_all,nbins,MPI_DOUBLE,MPI_SUM,world);
461 MPI_Allreduce(Pphi_temp,Pphi_all,nbins,MPI_DOUBLE,MPI_SUM,world);
462 MPI_Allreduce(Pz_temp,Pz_all,nbins,MPI_DOUBLE,MPI_SUM,world);
463
464 // populate array
465 for (ibin = 0; ibin < nbins; ibin++) {
466 array[ibin][0]=R[ibin];
467 array[ibin][2]=Pr_all[ibin]*nktv2p;
468 array[ibin][3]=Pphi_all[ibin]*nktv2p;
469 array[ibin][4]=Pz_all[ibin]*nktv2p;
470 }
471
472 }
473
474 /* ----------------------------------------------------------------------
475 memory usage of data
476 ------------------------------------------------------------------------- */
477
memory_usage()478 double ComputePressureCyl::memory_usage()
479 {
480 double bytes =
481 (3.0*(double)nphi + 16.0*(double)nbins+5.0*(double)nbins) * sizeof(double);
482 return bytes;
483 }
484