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 "ntopo_bond_partial.h"
16
17 #include "atom.h"
18 #include "force.h"
19 #include "domain.h"
20 #include "update.h"
21 #include "output.h"
22 #include "thermo.h"
23 #include "memory.h"
24 #include "error.h"
25
26
27 using namespace LAMMPS_NS;
28
29 #define DELTA 10000
30
31 /* ---------------------------------------------------------------------- */
32
NTopoBondPartial(LAMMPS * lmp)33 NTopoBondPartial::NTopoBondPartial(LAMMPS *lmp) : NTopo(lmp)
34 {
35 allocate_bond();
36 }
37
38 /* ---------------------------------------------------------------------- */
39
build()40 void NTopoBondPartial::build()
41 {
42 int i,m,atom1;
43
44 int nlocal = atom->nlocal;
45 int *num_bond = atom->num_bond;
46 tagint **bond_atom = atom->bond_atom;
47 int **bond_type = atom->bond_type;
48 tagint *tag = atom->tag;
49 int newton_bond = force->newton_bond;
50
51 int lostbond = output->thermo->lostbond;
52 int nmissing = 0;
53 nbondlist = 0;
54
55 for (i = 0; i < nlocal; i++)
56 for (m = 0; m < num_bond[i]; m++) {
57 if (bond_type[i][m] <= 0) continue;
58 atom1 = atom->map(bond_atom[i][m]);
59 if (atom1 == -1) {
60 nmissing++;
61 if (lostbond == Thermo::ERROR)
62 error->one(FLERR,"Bond atoms {} {} missing on "
63 "proc {} at step {}",tag[i],
64 bond_atom[i][m],me,update->ntimestep);
65 continue;
66 }
67 atom1 = domain->closest_image(i,atom1);
68 if (newton_bond || i < atom1) {
69 if (nbondlist == maxbond) {
70 maxbond += DELTA;
71 memory->grow(bondlist,maxbond,3,"neigh_topo:bondlist");
72 }
73 bondlist[nbondlist][0] = i;
74 bondlist[nbondlist][1] = atom1;
75 bondlist[nbondlist][2] = bond_type[i][m];
76 nbondlist++;
77 }
78 }
79
80 if (cluster_check) bond_check();
81 if (lostbond == Thermo::IGNORE) return;
82
83 int all;
84 MPI_Allreduce(&nmissing,&all,1,MPI_INT,MPI_SUM,world);
85 if (all && (me == 0))
86 error->warning(FLERR,"Bond atoms missing at step {}",update->ntimestep);
87 }
88