1 /* ----------------------------------------------------------------------
2     This is the
3 
4     ██╗     ██╗ ██████╗  ██████╗  ██████╗ ██╗  ██╗████████╗███████╗
5     ██║     ██║██╔════╝ ██╔════╝ ██╔════╝ ██║  ██║╚══██╔══╝██╔════╝
6     ██║     ██║██║  ███╗██║  ███╗██║  ███╗███████║   ██║   ███████╗
7     ██║     ██║██║   ██║██║   ██║██║   ██║██╔══██║   ██║   ╚════██║
8     ███████╗██║╚██████╔╝╚██████╔╝╚██████╔╝██║  ██║   ██║   ███████║
9     ╚══════╝╚═╝ ╚═════╝  ╚═════╝  ╚═════╝ ╚═╝  ╚═╝   ╚═╝   ╚══════╝®
10 
11     DEM simulation engine, released by
12     DCS Computing Gmbh, Linz, Austria
13     http://www.dcs-computing.com, office@dcs-computing.com
14 
15     LIGGGHTS® is part of CFDEM®project:
16     http://www.liggghts.com | http://www.cfdem.com
17 
18     Core developer and main author:
19     Christoph Kloss, christoph.kloss@dcs-computing.com
20 
21     LIGGGHTS® is open-source, distributed under the terms of the GNU Public
22     License, version 2 or later. It is distributed in the hope that it will
23     be useful, but WITHOUT ANY WARRANTY; without even the implied warranty
24     of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. You should have
25     received a copy of the GNU General Public License along with LIGGGHTS®.
26     If not, see http://www.gnu.org/licenses . See also top-level README
27     and LICENSE files.
28 
29     LIGGGHTS® and CFDEM® are registered trade marks of DCS Computing GmbH,
30     the producer of the LIGGGHTS® software and the CFDEM®coupling software
31     See http://www.cfdem.com/terms-trademark-policy for details.
32 
33 -------------------------------------------------------------------------
34     Contributing author and copyright for this file:
35     This file is from LAMMPS
36     LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
37     http://lammps.sandia.gov, Sandia National Laboratories
38     Steve Plimpton, sjplimp@sandia.gov
39 
40     Copyright (2003) Sandia Corporation.  Under the terms of Contract
41     DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
42     certain rights in this software.  This software is distributed under
43     the GNU General Public License.
44 ------------------------------------------------------------------------- */
45 
46 #include <cmath>
47 #include <string.h>
48 #include "ctype.h"
49 #include "angle_hybrid.h"
50 #include "atom.h"
51 #include "neighbor.h"
52 #include "domain.h"
53 #include "comm.h"
54 #include "force.h"
55 #include "memory.h"
56 #include "error.h"
57 
58 using namespace LAMMPS_NS;
59 
60 #define EXTRA 1000
61 
62 /* ---------------------------------------------------------------------- */
63 
AngleHybrid(LAMMPS * lmp)64 AngleHybrid::AngleHybrid(LAMMPS *lmp) : Angle(lmp)
65 {
66   writedata = 0;
67   nstyles = 0;
68 }
69 
70 /* ---------------------------------------------------------------------- */
71 
~AngleHybrid()72 AngleHybrid::~AngleHybrid()
73 {
74   if (nstyles) {
75     for (int i = 0; i < nstyles; i++) delete styles[i];
76     delete [] styles;
77     for (int i = 0; i < nstyles; i++) delete [] keywords[i];
78     delete [] keywords;
79   }
80 
81   if (allocated) {
82     memory->destroy(setflag);
83     memory->destroy(map);
84     delete [] nanglelist;
85     delete [] maxangle;
86     for (int i = 0; i < nstyles; i++)
87       memory->destroy(anglelist[i]);
88     delete [] anglelist;
89   }
90 }
91 
92 /* ---------------------------------------------------------------------- */
93 
compute(int eflag,int vflag)94 void AngleHybrid::compute(int eflag, int vflag)
95 {
96   int i,j,m,n;
97 
98   // save ptrs to original anglelist
99 
100   int nanglelist_orig = neighbor->nanglelist;
101   int **anglelist_orig = neighbor->anglelist;
102 
103   // if this is re-neighbor step, create sub-style anglelists
104   // nanglelist[] = length of each sub-style list
105   // realloc sub-style anglelist if necessary
106   // load sub-style anglelist with 4 values from original anglelist
107 
108   if (neighbor->ago == 0) {
109     for (m = 0; m < nstyles; m++) nanglelist[m] = 0;
110     for (i = 0; i < nanglelist_orig; i++) {
111       m = map[anglelist_orig[i][3]];
112       if (m >= 0) nanglelist[m]++;
113     }
114     for (m = 0; m < nstyles; m++) {
115       if (nanglelist[m] > maxangle[m]) {
116         memory->destroy(anglelist[m]);
117         maxangle[m] = nanglelist[m] + EXTRA;
118         memory->create(anglelist[m],maxangle[m],4,"angle_hybrid:anglelist");
119       }
120       nanglelist[m] = 0;
121     }
122     for (i = 0; i < nanglelist_orig; i++) {
123       m = map[anglelist_orig[i][3]];
124       if (m < 0) continue;
125       n = nanglelist[m];
126       anglelist[m][n][0] = anglelist_orig[i][0];
127       anglelist[m][n][1] = anglelist_orig[i][1];
128       anglelist[m][n][2] = anglelist_orig[i][2];
129       anglelist[m][n][3] = anglelist_orig[i][3];
130       nanglelist[m]++;
131     }
132   }
133 
134   // call each sub-style's compute function
135   // set neighbor->anglelist to sub-style anglelist before call
136   // accumulate sub-style global/peratom energy/virial in hybrid
137 
138   if (eflag || vflag) ev_setup(eflag,vflag);
139   else evflag = 0;
140 
141   for (m = 0; m < nstyles; m++) {
142     neighbor->nanglelist = nanglelist[m];
143     neighbor->anglelist = anglelist[m];
144 
145     styles[m]->compute(eflag,vflag);
146 
147     if (eflag_global) energy += styles[m]->energy;
148     if (vflag_global)
149       for (n = 0; n < 6; n++) virial[n] += styles[m]->virial[n];
150     if (eflag_atom) {
151       n = atom->nlocal;
152       if (force->newton_bond) n += atom->nghost;
153       double *eatom_substyle = styles[m]->eatom;
154       for (i = 0; i < n; i++) eatom[i] += eatom_substyle[i];
155     }
156     if (vflag_atom) {
157       n = atom->nlocal;
158       if (force->newton_bond) n += atom->nghost;
159       double **vatom_substyle = styles[m]->vatom;
160       for (i = 0; i < n; i++)
161         for (j = 0; j < 6; j++)
162           vatom[i][j] += vatom_substyle[i][j];
163     }
164   }
165 
166   // restore ptrs to original anglelist
167 
168   neighbor->nanglelist = nanglelist_orig;
169   neighbor->anglelist = anglelist_orig;
170 }
171 
172 /* ---------------------------------------------------------------------- */
173 
allocate()174 void AngleHybrid::allocate()
175 {
176   allocated = 1;
177   int n = atom->nangletypes;
178 
179   memory->create(map,n+1,"angle:map");
180   memory->create(setflag,n+1,"angle:setflag");
181   for (int i = 1; i <= n; i++) setflag[i] = 0;
182 
183   nanglelist = new int[nstyles];
184   maxangle = new int[nstyles];
185   anglelist = new int**[nstyles];
186   for (int m = 0; m < nstyles; m++) maxangle[m] = 0;
187   for (int m = 0; m < nstyles; m++) anglelist[m] = NULL;
188 }
189 
190 /* ----------------------------------------------------------------------
191    create one angle style for each arg in list
192 ------------------------------------------------------------------------- */
193 
settings(int narg,char ** arg)194 void AngleHybrid::settings(int narg, char **arg)
195 {
196   int i,m,istyle;
197 
198   if (narg < 1) error->all(FLERR,"Illegal angle_style command");
199 
200   // delete old lists, since cannot just change settings
201 
202   if (nstyles) {
203     for (int i = 0; i < nstyles; i++) delete styles[i];
204     delete [] styles;
205     for (int i = 0; i < nstyles; i++) delete [] keywords[i];
206     delete [] keywords;
207   }
208 
209   if (allocated) {
210     memory->destroy(setflag);
211     memory->destroy(map);
212     delete [] nanglelist;
213     delete [] maxangle;
214     for (int i = 0; i < nstyles; i++)
215       memory->destroy(anglelist[i]);
216     delete [] anglelist;
217   }
218   allocated = 0;
219 
220   // count sub-styles by skipping numeric args
221   // one exception is 1st arg of style "table", which is non-numeric word
222   // need a better way to skip these exceptions
223 
224   nstyles = 0;
225   i = 0;
226   while (i < narg) {
227     if (strcmp(arg[i],"table") == 0) i++;
228     i++;
229     while (i < narg && !isalpha(arg[i][0])) i++;
230     nstyles++;
231   }
232 
233   // allocate list of sub-styles
234 
235   styles = new Angle*[nstyles];
236   keywords = new char*[nstyles];
237 
238   // allocate each sub-style and call its settings() with subset of args
239   // define subset of args for a sub-style by skipping numeric args
240   // one exception is 1st arg of style "table", which is non-numeric
241   // need a better way to skip these exceptions
242 
243   int dummy;
244   nstyles = 0;
245   i = 0;
246 
247   while (i < narg) {
248     for (m = 0; m < nstyles; m++)
249       if (strcmp(arg[i],keywords[m]) == 0)
250         error->all(FLERR,"Angle style hybrid cannot use "
251                    "same angle style twice");
252     if (strcmp(arg[i],"hybrid") == 0)
253       error->all(FLERR,"Angle style hybrid cannot have hybrid as an argument");
254     if (strcmp(arg[i],"none") == 0)
255       error->all(FLERR,"Angle style hybrid cannot have none as an argument");
256     styles[nstyles] = force->new_angle(arg[i],lmp->suffix,dummy);
257     keywords[nstyles] = new char[strlen(arg[i])+1];
258     strcpy(keywords[nstyles],arg[i]);
259     istyle = i;
260     if (strcmp(arg[i],"table") == 0) i++;
261     i++;
262     while (i < narg && !isalpha(arg[i][0])) i++;
263     styles[nstyles]->settings(i-istyle-1,&arg[istyle+1]);
264     nstyles++;
265   }
266 }
267 
268 /* ----------------------------------------------------------------------
269    set coeffs for one type
270 ---------------------------------------------------------------------- */
271 
coeff(int narg,char ** arg)272 void AngleHybrid::coeff(int narg, char **arg)
273 {
274   if (!allocated) allocate();
275 
276   int ilo,ihi;
277   force->bounds(arg[0],atom->nangletypes,ilo,ihi);
278 
279   // 2nd arg = angle sub-style name
280   // allow for "none" or "skip" as valid sub-style name
281 
282   int m;
283   for (m = 0; m < nstyles; m++)
284     if (strcmp(arg[1],keywords[m]) == 0) break;
285 
286   int none = 0;
287   int skip = 0;
288   if (m == nstyles) {
289     if (strcmp(arg[1],"none") == 0) none = 1;
290     else if (strcmp(arg[1],"skip") == 0) none = skip = 1;
291     else error->all(FLERR,"Angle coeff for hybrid has invalid style");
292   }
293 
294   // move 1st arg to 2nd arg
295   // just copy ptrs, since arg[] points into original input line
296 
297   arg[1] = arg[0];
298 
299   // invoke sub-style coeff() starting with 1st arg
300 
301   if (!none) styles[m]->coeff(narg-1,&arg[1]);
302 
303   // set setflag and which type maps to which sub-style
304   // if sub-style is skip: auxiliary class2 setting in data file so ignore
305   // if sub-style is none: set hybrid setflag, wipe out map
306 
307   for (int i = ilo; i <= ihi; i++) {
308     if (skip) continue;
309     else if (none) {
310       setflag[i] = 1;
311       map[i] = -1;
312     } else {
313       setflag[i] = styles[m]->setflag[i];
314       map[i] = m;
315     }
316   }
317 }
318 
319 /* ----------------------------------------------------------------------
320    run angle style specific initialization
321 ------------------------------------------------------------------------- */
322 
init_style()323 void AngleHybrid::init_style()
324 {
325   for (int m = 0; m < nstyles; m++)
326     if (styles[m]) styles[m]->init_style();
327 }
328 
329 /* ----------------------------------------------------------------------
330    return an equilbrium angle length
331 ------------------------------------------------------------------------- */
332 
equilibrium_angle(int i)333 double AngleHybrid::equilibrium_angle(int i)
334 {
335   if (map[i] < 0)
336     error->one(FLERR,"Invoked angle equil angle on angle style none");
337   return styles[map[i]]->equilibrium_angle(i);
338 }
339 
340 /* ----------------------------------------------------------------------
341    proc 0 writes to restart file
342 ------------------------------------------------------------------------- */
343 
write_restart(FILE * fp)344 void AngleHybrid::write_restart(FILE *fp)
345 {
346   fwrite(&nstyles,sizeof(int),1,fp);
347 
348   int n;
349   for (int m = 0; m < nstyles; m++) {
350     n = strlen(keywords[m]) + 1;
351     fwrite(&n,sizeof(int),1,fp);
352     fwrite(keywords[m],sizeof(char),n,fp);
353   }
354 }
355 
356 /* ----------------------------------------------------------------------
357    proc 0 reads from restart file, bcasts
358 ------------------------------------------------------------------------- */
359 
read_restart(FILE * fp)360 void AngleHybrid::read_restart(FILE *fp)
361 {
362   int me = comm->me;
363   if (me == 0) fread(&nstyles,sizeof(int),1,fp);
364   MPI_Bcast(&nstyles,1,MPI_INT,0,world);
365   styles = new Angle*[nstyles];
366   keywords = new char*[nstyles];
367 
368   allocate();
369 
370   int n,dummy;
371   for (int m = 0; m < nstyles; m++) {
372     if (me == 0) fread(&n,sizeof(int),1,fp);
373     MPI_Bcast(&n,1,MPI_INT,0,world);
374     keywords[m] = new char[n];
375     if (me == 0) fread(keywords[m],sizeof(char),n,fp);
376     MPI_Bcast(keywords[m],n,MPI_CHAR,0,world);
377     styles[m] = force->new_angle(keywords[m],lmp->suffix,dummy);
378   }
379 }
380 
381 /* ---------------------------------------------------------------------- */
382 
single(int type,int i1,int i2,int i3)383 double AngleHybrid::single(int type, int i1, int i2, int i3)
384 {
385   if (map[type] < 0) error->one(FLERR,"Invoked angle single on angle style none");
386   return styles[map[type]]->single(type,i1,i2,i3);
387 }
388 
389 /* ----------------------------------------------------------------------
390    memory usage
391 ------------------------------------------------------------------------- */
392 
memory_usage()393 double AngleHybrid::memory_usage()
394 {
395   double bytes = maxeatom * sizeof(double);
396   bytes += maxvatom*6 * sizeof(double);
397   for (int m = 0; m < nstyles; m++) bytes += maxangle[m]*4 * sizeof(int);
398   for (int m = 0; m < nstyles; m++)
399     if (styles[m]) bytes += styles[m]->memory_usage();
400   return bytes;
401 }
402