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    This software is distributed under the GNU General Public License.
8 
9    See the README file in the top-level LAMMPS directory.
10 ------------------------------------------------------------------------- */
11 
12 /* ----------------------------------------------------------------------
13    Contributing author: Axel Kohlmeyer (Temple U)
14 ------------------------------------------------------------------------- */
15 
16 #include "pair_brownian_omp.h"
17 
18 #include "atom.h"
19 #include "comm.h"
20 #include "domain.h"
21 #include "fix_wall.h"
22 #include "force.h"
23 #include "input.h"
24 #include "math_const.h"
25 #include "math_special.h"
26 #include "neigh_list.h"
27 #include "random_mars.h"
28 #include "suffix.h"
29 #include "update.h"
30 #include "variable.h"
31 
32 #include <cmath>
33 
34 #include "omp_compat.h"
35 using namespace LAMMPS_NS;
36 using namespace MathConst;
37 using namespace MathSpecial;
38 
39 #define EPSILON 1.0e-10
40 
41 // same as fix_wall.cpp
42 
43 enum{EDGE,CONSTANT,VARIABLE};
44 
45 /* ---------------------------------------------------------------------- */
46 
PairBrownianOMP(LAMMPS * lmp)47 PairBrownianOMP::PairBrownianOMP(LAMMPS *lmp) :
48   PairBrownian(lmp), ThrOMP(lmp, THR_PAIR)
49 {
50   suffix_flag |= Suffix::OMP;
51   respa_enable = 0;
52   random_thr = nullptr;
53   nthreads = 0;
54 }
55 
56 /* ---------------------------------------------------------------------- */
57 
~PairBrownianOMP()58 PairBrownianOMP::~PairBrownianOMP()
59 {
60   if (random_thr) {
61     for (int i=1; i < nthreads; ++i)
62       delete random_thr[i];
63 
64     delete[] random_thr;
65     random_thr = nullptr;
66   }
67 }
68 
69 /* ---------------------------------------------------------------------- */
70 
compute(int eflag,int vflag)71 void PairBrownianOMP::compute(int eflag, int vflag)
72 {
73   ev_init(eflag,vflag);
74 
75   const int nall = atom->nlocal + atom->nghost;
76   const int inum = list->inum;
77 
78   // This section of code adjusts R0/RT0/RS0 if necessary due to changes
79   // in the volume fraction as a result of fix deform or moving walls
80 
81   double dims[3], wallcoord;
82   if (flagVF) // Flag for volume fraction corrections
83     if (flagdeform || flagwall == 2) { // Possible changes in volume fraction
84       if (flagdeform && !flagwall)
85         for (int j = 0; j < 3; j++)
86           dims[j] = domain->prd[j];
87       else if (flagwall == 2 || (flagdeform && flagwall == 1)) {
88         double wallhi[3], walllo[3];
89         for (int j = 0; j < 3; j++) {
90           wallhi[j] = domain->prd[j];
91           walllo[j] = 0;
92         }
93         for (int m = 0; m < wallfix->nwall; m++) {
94           int dim = wallfix->wallwhich[m] / 2;
95           int side = wallfix->wallwhich[m] % 2;
96           if (wallfix->xstyle[m] == VARIABLE) {
97             wallcoord = input->variable->compute_equal(wallfix->xindex[m]);
98           }
99           else wallcoord = wallfix->coord0[m];
100           if (side == 0) walllo[dim] = wallcoord;
101           else wallhi[dim] = wallcoord;
102         }
103         for (int j = 0; j < 3; j++)
104           dims[j] = wallhi[j] - walllo[j];
105       }
106       double vol_T = dims[0]*dims[1]*dims[2];
107       double vol_f = vol_P/vol_T;
108       if (flaglog == 0) {
109         R0  = 6*MY_PI*mu*rad*(1.0 + 2.16*vol_f);
110         RT0 = 8*MY_PI*mu*cube(rad);
111         //RS0 = 20.0/3.0*MY_PI*mu*pow(rad,3)*(1.0 + 3.33*vol_f + 2.80*vol_f*vol_f);
112       } else {
113         R0  = 6*MY_PI*mu*rad*(1.0 + 2.725*vol_f - 6.583*vol_f*vol_f);
114         RT0 = 8*MY_PI*mu*cube(rad)*(1.0 + 0.749*vol_f - 2.469*vol_f*vol_f);
115         //RS0 = 20.0/3.0*MY_PI*mu*pow(rad,3)*(1.0 + 3.64*vol_f - 6.95*vol_f*vol_f);
116       }
117     }
118 
119   // number of threads has changed. reallocate pool of pRNGs
120   if (nthreads != comm->nthreads) {
121     if (random_thr) {
122       for (int i=1; i < nthreads; ++i)
123         delete random_thr[i];
124 
125       delete[] random_thr;
126     }
127 
128     nthreads = comm->nthreads;
129     random_thr = new RanMars*[nthreads];
130     for (int i=1; i < nthreads; ++i)
131       random_thr[i] = nullptr;
132 
133     // to ensure full compatibility with the serial Brownian style
134     // we use is random number generator instance for thread 0
135     random_thr[0] = random;
136   }
137 
138 #if defined(_OPENMP)
139 #pragma omp parallel LMP_DEFAULT_NONE LMP_SHARED(eflag,vflag)
140 #endif
141   {
142     int ifrom, ito, tid;
143 
144     loop_setup_thr(ifrom, ito, tid, inum, nthreads);
145     ThrData *thr = fix->get_thr(tid);
146     thr->timer(Timer::START);
147     ev_setup_thr(eflag, vflag, nall, eatom, vatom, nullptr, thr);
148 
149     // generate a random number generator instance for
150     // all threads != 0. make sure we use unique seeds.
151     if ((tid > 0) && (random_thr[tid] == nullptr))
152       random_thr[tid] = new RanMars(Pair::lmp, seed + comm->me
153                                     + comm->nprocs*tid);
154 
155     if (flaglog) {
156       if (evflag) {
157         if (force->newton_pair) eval<1,1,1>(ifrom, ito, thr);
158         else eval<1,1,0>(ifrom, ito, thr);
159       } else {
160         if (force->newton_pair) eval<1,0,1>(ifrom, ito, thr);
161         else eval<1,0,0>(ifrom, ito, thr);
162       }
163     } else {
164       if (evflag) {
165         if (force->newton_pair) eval<0,1,1>(ifrom, ito, thr);
166         else eval<0,1,0>(ifrom, ito, thr);
167       } else {
168         if (force->newton_pair) eval<0,0,1>(ifrom, ito, thr);
169         else eval<0,0,0>(ifrom, ito, thr);
170       }
171     }
172 
173     thr->timer(Timer::PAIR);
174     reduce_thr(this, eflag, vflag, thr);
175   } // end of omp parallel region
176 }
177 
178 template <int FLAGLOG, int EVFLAG, int NEWTON_PAIR>
eval(int iifrom,int iito,ThrData * const thr)179 void PairBrownianOMP::eval(int iifrom, int iito, ThrData * const thr)
180 {
181   int i,j,ii,jj,jnum,itype,jtype;
182   double xtmp,ytmp,ztmp,delx,dely,delz,fx,fy,fz,tx,ty,tz;
183   double rsq,r,h_sep,radi;
184   int *ilist,*jlist,*numneigh,**firstneigh;
185 
186   const double * const * const x = atom->x;
187   double * const * const f = thr->get_f();
188   double * const * const torque = thr->get_torque();
189   const double * const radius = atom->radius;
190   const int * const type = atom->type;
191   const int nlocal = atom->nlocal;
192 
193   RanMars &rng = *random_thr[thr->get_tid()];
194 
195   double vxmu2f = force->vxmu2f;
196   double randr;
197   double prethermostat;
198   double xl[3],a_sq,a_sh,a_pu,Fbmag;
199   double p1[3],p2[3],p3[3];
200   int overlaps = 0;
201 
202   // scale factor for Brownian moments
203 
204   prethermostat = sqrt(24.0*force->boltz*t_target/update->dt);
205   prethermostat *= sqrt(force->vxmu2f/force->ftm2v/force->mvv2e);
206 
207   ilist = list->ilist;
208   numneigh = list->numneigh;
209   firstneigh = list->firstneigh;
210 
211   // loop over neighbors of my atoms
212 
213   for (ii = iifrom; ii < iito; ++ii) {
214     i = ilist[ii];
215     xtmp = x[i][0];
216     ytmp = x[i][1];
217     ztmp = x[i][2];
218     itype = type[i];
219     radi = radius[i];
220     jlist = firstneigh[i];
221     jnum = numneigh[i];
222 
223     // FLD contribution to force and torque due to isotropic terms
224 
225     if (flagfld) {
226       f[i][0] += prethermostat*sqrt(R0)*(rng.uniform()-0.5);
227       f[i][1] += prethermostat*sqrt(R0)*(rng.uniform()-0.5);
228       f[i][2] += prethermostat*sqrt(R0)*(rng.uniform()-0.5);
229       if (FLAGLOG) {
230         torque[i][0] += prethermostat*sqrt(RT0)*(rng.uniform()-0.5);
231         torque[i][1] += prethermostat*sqrt(RT0)*(rng.uniform()-0.5);
232         torque[i][2] += prethermostat*sqrt(RT0)*(rng.uniform()-0.5);
233       }
234     }
235 
236     if (!flagHI) continue;
237 
238     for (jj = 0; jj < jnum; jj++) {
239       j = jlist[jj];
240       j &= NEIGHMASK;
241 
242       delx = xtmp - x[j][0];
243       dely = ytmp - x[j][1];
244       delz = ztmp - x[j][2];
245       rsq = delx*delx + dely*dely + delz*delz;
246       jtype = type[j];
247 
248       if (rsq < cutsq[itype][jtype]) {
249         r = sqrt(rsq);
250 
251         // scalar resistances a_sq and a_sh
252 
253         h_sep = r - 2.0*radi;
254 
255         // check for overlaps
256 
257         if (h_sep < 0.0) overlaps++;
258 
259         // if less than minimum gap, use minimum gap instead
260 
261         if (r < cut_inner[itype][jtype])
262           h_sep = cut_inner[itype][jtype] - 2.0*radi;
263 
264         // scale h_sep by radi
265 
266         h_sep = h_sep/radi;
267 
268         // scalar resistances
269 
270         if (FLAGLOG) {
271           a_sq = 6.0*MY_PI*mu*radi*(1.0/4.0/h_sep + 9.0/40.0*log(1.0/h_sep));
272           a_sh = 6.0*MY_PI*mu*radi*(1.0/6.0*log(1.0/h_sep));
273           a_pu = 8.0*MY_PI*mu*cube(radi)*(3.0/160.0*log(1.0/h_sep));
274         } else
275           a_sq = 6.0*MY_PI*mu*radi*(1.0/4.0/h_sep);
276 
277         // generate the Pairwise Brownian Force: a_sq
278 
279         Fbmag = prethermostat*sqrt(a_sq);
280 
281         // generate a random number
282 
283         randr = rng.uniform()-0.5;
284 
285         // contribution due to Brownian motion
286 
287         fx = Fbmag*randr*delx/r;
288         fy = Fbmag*randr*dely/r;
289         fz = Fbmag*randr*delz/r;
290 
291         // add terms due to a_sh
292 
293         if (FLAGLOG) {
294 
295           // generate two orthogonal vectors to the line of centers
296 
297           p1[0] = delx/r; p1[1] = dely/r; p1[2] = delz/r;
298           set_3_orthogonal_vectors(p1,p2,p3);
299 
300           // magnitude
301 
302           Fbmag = prethermostat*sqrt(a_sh);
303 
304           // force in each of the two directions
305 
306           randr = rng.uniform()-0.5;
307           fx += Fbmag*randr*p2[0];
308           fy += Fbmag*randr*p2[1];
309           fz += Fbmag*randr*p2[2];
310 
311           randr = rng.uniform()-0.5;
312           fx += Fbmag*randr*p3[0];
313           fy += Fbmag*randr*p3[1];
314           fz += Fbmag*randr*p3[2];
315         }
316 
317         // scale forces to appropriate units
318 
319         fx = vxmu2f*fx;
320         fy = vxmu2f*fy;
321         fz = vxmu2f*fz;
322 
323         // sum to total force
324 
325         f[i][0] -= fx;
326         f[i][1] -= fy;
327         f[i][2] -= fz;
328 
329         if (NEWTON_PAIR || j < nlocal) {
330           //randr = rng.uniform()-0.5;
331           //fx = Fbmag*randr*delx/r;
332           //fy = Fbmag*randr*dely/r;
333           //fz = Fbmag*randr*delz/r;
334 
335           f[j][0] += fx;
336           f[j][1] += fy;
337           f[j][2] += fz;
338         }
339 
340         // torque due to the Brownian Force
341 
342         if (FLAGLOG) {
343 
344           // location of the point of closest approach on I from its center
345 
346           xl[0] = -delx/r*radi;
347           xl[1] = -dely/r*radi;
348           xl[2] = -delz/r*radi;
349 
350           // torque = xl_cross_F
351 
352           tx = xl[1]*fz - xl[2]*fy;
353           ty = xl[2]*fx - xl[0]*fz;
354           tz = xl[0]*fy - xl[1]*fx;
355 
356           // torque is same on both particles
357 
358           torque[i][0] -= tx;
359           torque[i][1] -= ty;
360           torque[i][2] -= tz;
361 
362           if (NEWTON_PAIR || j < nlocal) {
363             torque[j][0] -= tx;
364             torque[j][1] -= ty;
365             torque[j][2] -= tz;
366           }
367 
368           // torque due to a_pu
369 
370           Fbmag = prethermostat*sqrt(a_pu);
371 
372           // force in each direction
373 
374           randr = rng.uniform()-0.5;
375           tx = Fbmag*randr*p2[0];
376           ty = Fbmag*randr*p2[1];
377           tz = Fbmag*randr*p2[2];
378 
379           randr = rng.uniform()-0.5;
380           tx += Fbmag*randr*p3[0];
381           ty += Fbmag*randr*p3[1];
382           tz += Fbmag*randr*p3[2];
383 
384           // torque has opposite sign on two particles
385 
386           torque[i][0] -= tx;
387           torque[i][1] -= ty;
388           torque[i][2] -= tz;
389 
390           if (NEWTON_PAIR || j < nlocal) {
391             torque[j][0] += tx;
392             torque[j][1] += ty;
393             torque[j][2] += tz;
394           }
395         }
396 
397         if (EVFLAG) ev_tally_xyz(i,j,nlocal,NEWTON_PAIR,
398                                  0.0,0.0,-fx,-fy,-fz,delx,dely,delz);
399       }
400     }
401   }
402 }
403 
404 /* ---------------------------------------------------------------------- */
405 
memory_usage()406 double PairBrownianOMP::memory_usage()
407 {
408   double bytes = memory_usage_thr();
409   bytes += PairBrownian::memory_usage();
410   bytes += (double)nthreads * sizeof(RanMars*);
411   bytes += (double)nthreads * sizeof(RanMars);
412 
413   return bytes;
414 }
415