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_dihedral_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 
NTopoDihedralPartial(LAMMPS * lmp)33 NTopoDihedralPartial::NTopoDihedralPartial(LAMMPS *lmp) :
34   NTopo(lmp)
35 {
36   allocate_dihedral();
37 }
38 
39 /* ---------------------------------------------------------------------- */
40 
build()41 void NTopoDihedralPartial::build()
42 {
43   int i,m,atom1,atom2,atom3,atom4;
44 
45   int nlocal = atom->nlocal;
46   int *num_dihedral = atom->num_dihedral;
47   tagint **dihedral_atom1 = atom->dihedral_atom1;
48   tagint **dihedral_atom2 = atom->dihedral_atom2;
49   tagint **dihedral_atom3 = atom->dihedral_atom3;
50   tagint **dihedral_atom4 = atom->dihedral_atom4;
51   int **dihedral_type = atom->dihedral_type;
52   int newton_bond = force->newton_bond;
53 
54   int lostbond = output->thermo->lostbond;
55   int nmissing = 0;
56   ndihedrallist = 0;
57 
58   for (i = 0; i < nlocal; i++)
59     for (m = 0; m < num_dihedral[i]; m++) {
60       if (dihedral_type[i][m] <= 0) continue;
61       atom1 = atom->map(dihedral_atom1[i][m]);
62       atom2 = atom->map(dihedral_atom2[i][m]);
63       atom3 = atom->map(dihedral_atom3[i][m]);
64       atom4 = atom->map(dihedral_atom4[i][m]);
65       if (atom1 == -1 || atom2 == -1 || atom3 == -1 || atom4 == -1) {
66         nmissing++;
67         if (lostbond == Thermo::ERROR)
68           error->one(FLERR,"Dihedral atoms {} {} {} {} missing on "
69                                        "proc {} at step {}",
70                                        dihedral_atom1[i][m],dihedral_atom2[i][m],
71                                        dihedral_atom3[i][m],dihedral_atom4[i][m],
72                                        me,update->ntimestep);
73         continue;
74       }
75       atom1 = domain->closest_image(i,atom1);
76       atom2 = domain->closest_image(i,atom2);
77       atom3 = domain->closest_image(i,atom3);
78       atom4 = domain->closest_image(i,atom4);
79       if (newton_bond ||
80           (i <= atom1 && i <= atom2 && i <= atom3 && i <= atom4)) {
81         if (ndihedrallist == maxdihedral) {
82           maxdihedral += DELTA;
83           memory->grow(dihedrallist,maxdihedral,5,"neigh_topo:dihedrallist");
84         }
85         dihedrallist[ndihedrallist][0] = atom1;
86         dihedrallist[ndihedrallist][1] = atom2;
87         dihedrallist[ndihedrallist][2] = atom3;
88         dihedrallist[ndihedrallist][3] = atom4;
89         dihedrallist[ndihedrallist][4] = dihedral_type[i][m];
90         ndihedrallist++;
91       }
92     }
93 
94   if (cluster_check) dihedral_check(ndihedrallist,dihedrallist);
95   if (lostbond == Thermo::IGNORE) return;
96 
97   int all;
98   MPI_Allreduce(&nmissing,&all,1,MPI_INT,MPI_SUM,world);
99   if (all && (me == 0))
100     error->warning(FLERR,"Dihedral atoms missing at step {}",update->ntimestep);
101 }
102