lammps-sys 0.6.0

Generates bindings to LAMMPS' C interface (with optional builds from source)
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
/* ----------------------------------------------------------------------
   LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
   http://lammps.sandia.gov, Sandia National Laboratories
   Steve Plimpton, sjplimp@sandia.gov

   Copyright (2003) Sandia Corporation.  Under the terms of Contract
   DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
   certain rights in this software.  This software is distributed under
   the GNU General Public License.

   See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */

/* ----------------------------------------------------------------------
   Contributing authors: Mike Salerno (NRL) added single methods
------------------------------------------------------------------------- */

#include "create_bonds.h"
#include <mpi.h>
#include <cstring>
#include "atom.h"
#include "domain.h"
#include "force.h"
#include "neighbor.h"
#include "neigh_request.h"
#include "neigh_list.h"
#include "comm.h"
#include "group.h"
#include "special.h"
#include "error.h"

using namespace LAMMPS_NS;

enum{MANY,SBOND,SANGLE,SDIHEDRAL};

/* ---------------------------------------------------------------------- */

CreateBonds::CreateBonds(LAMMPS *lmp) : Pointers(lmp) {}

/* ---------------------------------------------------------------------- */

void CreateBonds::command(int narg, char **arg)
{
  if (domain->box_exist == 0)
    error->all(FLERR,"Create_bonds command before simulation box is defined");
  if (atom->tag_enable == 0)
    error->all(FLERR,"Cannot use create_bonds unless atoms have IDs");
  if (atom->molecular != 1)
    error->all(FLERR,"Cannot use create_bonds with non-molecular system");

  if (narg < 4) error->all(FLERR,"Illegal create_bonds command");

  // parse args

  int style;

  int iarg = 0;
  if (strcmp(arg[0],"many") == 0) {
    style = MANY;
    if (narg != 6) error->all(FLERR,"Illegal create_bonds command");
    igroup = group->find(arg[1]);
    if (igroup == -1) error->all(FLERR,"Cannot find create_bonds group ID");
    group1bit = group->bitmask[igroup];
    igroup = group->find(arg[2]);
    if (igroup == -1) error->all(FLERR,"Cannot find create_bonds group ID");
    group2bit = group->bitmask[igroup];
    btype = force->inumeric(FLERR,arg[3]);
    rmin = force->numeric(FLERR,arg[4]);
    rmax = force->numeric(FLERR,arg[5]);
    if (rmin > rmax) error->all(FLERR,"Illegal create_bonds command");
    iarg = 6;
  } else if (strcmp(arg[0],"single/bond") == 0) {
    style = SBOND;
    if (narg < 4) error->all(FLERR,"Illegal create_bonds command");
    btype = force->inumeric(FLERR,arg[1]);
    batom1 = force->tnumeric(FLERR,arg[2]);
    batom2 = force->tnumeric(FLERR,arg[3]);
    if (batom1 == batom2)
      error->all(FLERR,"Illegal create_bonds command");
    iarg = 4;
  } else if (strcmp(arg[0],"single/angle") == 0) {
    style = SANGLE;
    if (narg < 5) error->all(FLERR,"Illegal create_bonds command");
    atype = force->inumeric(FLERR,arg[1]);
    aatom1 = force->tnumeric(FLERR,arg[2]);
    aatom2 = force->tnumeric(FLERR,arg[3]);
    aatom3 = force->tnumeric(FLERR,arg[4]);
    if ((aatom1 == aatom2) || (aatom1 == aatom3) || (aatom2 == aatom3))
      error->all(FLERR,"Illegal create_bonds command");
    iarg = 5;
  } else if (strcmp(arg[0],"single/dihedral") == 0) {
    style = SDIHEDRAL;
    if (narg < 6) error->all(FLERR,"Illegal create_bonds command");
    dtype = force->inumeric(FLERR,arg[1]);
    datom1 = force->tnumeric(FLERR,arg[2]);
    datom2 = force->tnumeric(FLERR,arg[3]);
    datom3 = force->tnumeric(FLERR,arg[4]);
    datom4 = force->tnumeric(FLERR,arg[5]);
    if ((datom1 == datom2) || (datom1 == datom3) || (datom1 == datom4) ||
        (datom2 == datom3) || (datom2 == datom4) || (datom3 == datom4))
      error->all(FLERR,"Illegal create_bonds command");
    iarg = 6;
  } else error->all(FLERR,"Illegal create_bonds command");

  // optional args

  int specialflag = 1;

  while (iarg < narg) {
    if (strcmp(arg[iarg],"special") == 0) {
      if (iarg+2 > narg) error->all(FLERR,"Illegal create_bonds command");
      if (strcmp(arg[iarg+1],"yes") == 0) specialflag = 1;
      else if (strcmp(arg[iarg+1],"no") == 0) specialflag = 0;
      else error->all(FLERR,"Illegal create_bonds command");
      iarg += 2;
    } else error->all(FLERR,"Illegal create_bonds command");
  }

  // error checks

  if (style == MANY) {
    if (btype <= 0 || btype > atom->nbondtypes)
      error->all(FLERR,"Invalid bond type in create_bonds command");
    if (specialflag == 0)
      error->all(FLERR,"Cannot use special no with create_bonds many");
  } else if (style == SBOND) {
    if (btype <= 0 || btype > atom->nbondtypes)
      error->all(FLERR,"Invalid bond type in create_bonds command");
  } else if (style == SANGLE) {
    if (atype <= 0 || atype > atom->nangletypes)
      error->all(FLERR,"Invalid angle type in create_bonds command");
  } else if (style == SDIHEDRAL) {
    if (dtype <= 0 || dtype > atom->ndihedraltypes)
      error->all(FLERR,"Invalid dihedral type in create_bonds command");
  }

  // invoke creation method

  if (style == MANY) many();
  else if (style == SBOND) single_bond();
  else if (style == SANGLE) single_angle();
  else if (style == SDIHEDRAL) single_dihedral();

  // trigger special list build

  if (specialflag) {
    Special special(lmp);
    special.build();
  }
}

/* ---------------------------------------------------------------------- */

void CreateBonds::many()
{
  double rminsq = rmin*rmin;
  double rmaxsq = rmax*rmax;

  // store state before bond creation

  bigint nbonds_previous = atom->nbonds;

  // request a full neighbor list for use by this command

  int irequest = neighbor->request(this);
  neighbor->requests[irequest]->pair = 0;
  neighbor->requests[irequest]->command = 1;
  neighbor->requests[irequest]->half = 0;
  neighbor->requests[irequest]->full = 1;
  neighbor->requests[irequest]->occasional = 1;
  neighbor->requests[irequest]->command_style = "create_bonds";

  // init entire system since comm->borders and neighbor->build is done
  // comm::init needs neighbor::init needs pair::init needs kspace::init, etc

  lmp->init();

  // error check on cutoff
  // if no pair style, neighbor list will be empty

  if (force->pair == NULL)
    error->all(FLERR,"Create_bonds requires a pair style be defined");
  if (rmax > neighbor->cutneighmax)
    error->all(FLERR,"Create_bonds max distance > neighbor cutoff");
  if (rmax > neighbor->cutneighmin && comm->me == 0)
    error->warning(FLERR,"Create_bonds max distance > minimum neighbor cutoff");

  // require special_bonds 1-2 weights = 0.0 and KSpace = NULL
  // so that already bonded atom pairs do not appear in neighbor list
  // otherwise with newton_bond = 1,
  //   would be hard to check if I-J bond already existed
  // note that with KSpace, pair with weight = 0 could still be in neigh list

  if (force->special_lj[1] != 0.0 || force->special_coul[1] != 0.0)
    error->all(FLERR,"Create_bonds command requires "
               "special_bonds 1-2 weights be 0.0");
  if (force->kspace)
    error->all(FLERR,"Create_bonds command requires "
               "no kspace_style be defined");

  // setup domain, communication and neighboring
  // acquire ghosts and build standard neighbor lists

  if (domain->triclinic) domain->x2lamda(atom->nlocal);
  domain->pbc();
  domain->reset_box();
  comm->setup();
  if (neighbor->style) neighbor->setup_bins();
  comm->exchange();
  comm->borders();
  if (domain->triclinic) domain->lamda2x(atom->nlocal+atom->nghost);
  neighbor->build(1);

  // build neighbor list this command needs based on earlier request

  NeighList *list = neighbor->lists[irequest];
  neighbor->build_one(list);

  // loop over all neighs of each atom
  // compute distance between two atoms consistently on both procs
  // add bond if group and distance criteria are met
  // check that bond list does not overflow

  tagint *tag = atom->tag;
  int *mask = atom->mask;
  double **x = atom->x;
  int *num_bond = atom->num_bond;
  int **bond_type = atom->bond_type;
  tagint **bond_atom = atom->bond_atom;
  double newton_bond = force->newton_bond;
  int nlocal = atom->nlocal;

  int i,j,ii,jj,inum,jnum,flag;
  double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
  int *ilist,*jlist,*numneigh,**firstneigh;

  inum = list->inum;
  ilist = list->ilist;
  numneigh = list->numneigh;
  firstneigh = list->firstneigh;

  for (ii = 0; ii < inum; ii++) {
    i = ilist[ii];
    xtmp = x[i][0];
    ytmp = x[i][1];
    ztmp = x[i][2];
    jlist = firstneigh[i];
    jnum = numneigh[i];

    for (jj = 0; jj < jnum; jj++) {
      j = jlist[jj];
      j &= NEIGHMASK;

      // only consider bond creation if I,J distance between 2 cutoffs
      // compute rsq identically on both I,J loop iterations
      // if I,J tags equal, do not bond atom to itself

      if (tag[i] < tag[j]) {
        delx = xtmp - x[j][0];
        dely = ytmp - x[j][1];
        delz = ztmp - x[j][2];
      } else if (tag[i] > tag[j]) {
        delx = x[j][0] - xtmp;
        dely = x[j][1] - ytmp;
        delz = x[j][2] - ztmp;
      } else continue;
      rsq = delx*delx + dely*dely + delz*delz;
      if (rsq < rminsq || rsq > rmaxsq) continue;

      // only consider bond creation if igroup and jgroup match I,J atoms

      flag = 0;
      if ((mask[i] & group1bit) && (mask[j] & group2bit)) flag = 1;
      if ((mask[i] & group2bit) && (mask[j] & group1bit)) flag = 1;
      if (!flag) continue;

      // create bond, check for overflow
      // on I,J loop iterations, store with 1 or 2 atoms based on newton_bond

      if (!newton_bond || tag[i] < tag[j]) {
        if (num_bond[i] == atom->bond_per_atom)
          error->one(FLERR,
                     "New bond exceeded bonds per atom in create_bonds");
        bond_type[i][num_bond[i]] = btype;
        bond_atom[i][num_bond[i]] = tag[j];
        num_bond[i]++;
      }
    }
  }

  // recount bonds

  bigint nbonds = 0;
  for (int i = 0; i < nlocal; i++) nbonds += num_bond[i];

  MPI_Allreduce(&nbonds,&atom->nbonds,1,MPI_LMP_BIGINT,MPI_SUM,world);
  if (!force->newton_bond) atom->nbonds /= 2;

  // print new bond count

  bigint nadd_bonds = atom->nbonds - nbonds_previous;

  if (comm->me == 0) {
    if (screen) {
      fprintf(screen,"Added " BIGINT_FORMAT
              " bonds, new total = " BIGINT_FORMAT "\n",
              nadd_bonds,atom->nbonds);
    }

    if (logfile) {
      fprintf(logfile,"Added " BIGINT_FORMAT
              " bonds, new total = " BIGINT_FORMAT "\n",
              nadd_bonds,atom->nbonds);
    }
  }
}

/* ---------------------------------------------------------------------- */

void CreateBonds::single_bond()
{
  int m;

  // check that 2 atoms exist

  const int nlocal = atom->nlocal;
  const int idx1 = atom->map(batom1);
  const int idx2 = atom->map(batom2);

  int count = 0;
  if ((idx1 >= 0) && (idx1 < nlocal)) count++;
  if ((idx2 >= 0) && (idx2 < nlocal)) count++;

  int allcount;
  MPI_Allreduce(&count,&allcount,1,MPI_INT,MPI_SUM,world);
  if (allcount != 2)
    error->all(FLERR,"Create_bonds single/bond atoms do not exist");

  // create bond once or 2x if newton_bond set

  int *num_bond = atom->num_bond;
  int **bond_type = atom->bond_type;
  tagint **bond_atom = atom->bond_atom;

  if ((m = idx1) >= 0) {
    if (num_bond[m] == atom->bond_per_atom)
      error->one(FLERR,"New bond exceeded bonds per atom in create_bonds");
    bond_type[m][num_bond[m]] = btype;
    bond_atom[m][num_bond[m]] = batom2;
    num_bond[m]++;
  }
  atom->nbonds++;

  if (force->newton_bond) return;

  if ((m = idx2) >= 0) {
    if (num_bond[m] == atom->bond_per_atom)
      error->one(FLERR,"New bond exceeded bonds per atom in create_bonds");
    bond_type[m][num_bond[m]] = btype;
    bond_atom[m][num_bond[m]] = batom1;
    num_bond[m]++;
  }
}

/* ---------------------------------------------------------------------- */

void CreateBonds::single_angle()
{
  int m;

  // check that 3 atoms exist

  const int nlocal = atom->nlocal;
  const int idx1 = atom->map(aatom1);
  const int idx2 = atom->map(aatom2);
  const int idx3 = atom->map(aatom3);

  int count = 0;
  if ((idx1 >= 0) && (idx1 < nlocal)) count++;
  if ((idx2 >= 0) && (idx2 < nlocal)) count++;
  if ((idx3 >= 0) && (idx3 < nlocal)) count++;

  int allcount;
  MPI_Allreduce(&count,&allcount,1,MPI_INT,MPI_SUM,world);
  if (allcount != 3)
    error->all(FLERR,"Create_bonds single/angle atoms do not exist");

  // create angle once or 3x if newton_bond set

  int *num_angle = atom->num_angle;
  int **angle_type = atom->angle_type;
  tagint **angle_atom1 = atom->angle_atom1;
  tagint **angle_atom2 = atom->angle_atom2;
  tagint **angle_atom3 = atom->angle_atom3;

  if ((m = idx2) >= 0) {
    if (num_angle[m] == atom->angle_per_atom)
      error->one(FLERR,"New angle exceeded angles per atom in create_bonds");
    angle_type[m][num_angle[m]] = atype;
    angle_atom1[m][num_angle[m]] = aatom1;
    angle_atom2[m][num_angle[m]] = aatom2;
    angle_atom3[m][num_angle[m]] = aatom3;
    num_angle[m]++;
  }
  atom->nangles++;

  if (force->newton_bond) return;

  if ((m = idx1) >= 0) {
    if (num_angle[m] == atom->angle_per_atom)
      error->one(FLERR,"New angle exceeded angles per atom in create_bonds");
    angle_type[m][num_angle[m]] = atype;
    angle_atom1[m][num_angle[m]] = aatom1;
    angle_atom2[m][num_angle[m]] = aatom2;
    angle_atom3[m][num_angle[m]] = aatom3;
    num_angle[m]++;
  }

  if ((m = idx3) >= 0) {
    if (num_angle[m] == atom->angle_per_atom)
      error->one(FLERR,"New angle exceeded angles per atom in create_bonds");
    angle_type[m][num_angle[m]] = atype;
    angle_atom1[m][num_angle[m]] = aatom1;
    angle_atom2[m][num_angle[m]] = aatom2;
    angle_atom3[m][num_angle[m]] = aatom3;
    num_angle[m]++;
  }
}

/* ---------------------------------------------------------------------- */

void CreateBonds::single_dihedral()
{
  int m;

  // check that 4 atoms exist

  const int nlocal = atom->nlocal;
  const int idx1 = atom->map(datom1);
  const int idx2 = atom->map(datom2);
  const int idx3 = atom->map(datom3);
  const int idx4 = atom->map(datom4);

  int count = 0;
  if ((idx1 >= 0) && (idx1 < nlocal)) count++;
  if ((idx2 >= 0) && (idx2 < nlocal)) count++;
  if ((idx3 >= 0) && (idx3 < nlocal)) count++;
  if ((idx4 >= 0) && (idx4 < nlocal)) count++;

  int allcount;
  MPI_Allreduce(&count,&allcount,1,MPI_INT,MPI_SUM,world);
  if (allcount != 4)
    error->all(FLERR,"Create_bonds single/dihedral atoms do not exist");

  // create bond once or 4x if newton_bond set

  int *num_dihedral = atom->num_dihedral;
  int **dihedral_type = atom->dihedral_type;
  tagint **dihedral_atom1 = atom->dihedral_atom1;
  tagint **dihedral_atom2 = atom->dihedral_atom2;
  tagint **dihedral_atom3 = atom->dihedral_atom3;
  tagint **dihedral_atom4 = atom->dihedral_atom4;

  if ((m = idx2) >= 0) {
    if (num_dihedral[m] == atom->dihedral_per_atom)
      error->one(FLERR,
                 "New dihedral exceeded dihedrals per atom in create_bonds");
    dihedral_type[m][num_dihedral[m]] = dtype;
    dihedral_atom1[m][num_dihedral[m]] = datom1;
    dihedral_atom2[m][num_dihedral[m]] = datom2;
    dihedral_atom3[m][num_dihedral[m]] = datom3;
    dihedral_atom4[m][num_dihedral[m]] = datom4;
    num_dihedral[m]++;
  }
  atom->ndihedrals++;

  if (force->newton_bond) return;

  if ((m = idx1) >= 0) {
    if (num_dihedral[m] == atom->dihedral_per_atom)
      error->one(FLERR,
                 "New dihedral exceeded dihedrals per atom in create_bonds");
    dihedral_type[m][num_dihedral[m]] = dtype;
    dihedral_atom1[m][num_dihedral[m]] = datom1;
    dihedral_atom2[m][num_dihedral[m]] = datom2;
    dihedral_atom3[m][num_dihedral[m]] = datom3;
    dihedral_atom4[m][num_dihedral[m]] = datom4;
    num_dihedral[m]++;
  }

  if ((m = idx3) >= 0) {
    if (num_dihedral[m] == atom->dihedral_per_atom)
      error->one(FLERR,
                 "New dihedral exceeded dihedrals per atom in create_bonds");
    dihedral_type[m][num_dihedral[m]] = dtype;
    dihedral_atom1[m][num_dihedral[m]] = datom1;
    dihedral_atom2[m][num_dihedral[m]] = datom2;
    dihedral_atom3[m][num_dihedral[m]] = datom3;
    dihedral_atom4[m][num_dihedral[m]] = datom4;
    num_dihedral[m]++;
  }

  if ((m = idx4) >= 0) {
    if (num_dihedral[m] == atom->dihedral_per_atom)
      error->one(FLERR,
                 "New dihedral exceeded dihedrals per atom in create_bonds");
    dihedral_type[m][num_dihedral[m]] = dtype;
    dihedral_atom1[m][num_dihedral[m]] = datom1;
    dihedral_atom2[m][num_dihedral[m]] = datom2;
    dihedral_atom3[m][num_dihedral[m]] = datom3;
    dihedral_atom4[m][num_dihedral[m]] = datom4;
    num_dihedral[m]++;
  }
}