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
|
/* ----------------------------------------------------------------------
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 author: Wan Liang (Chinese Academy of Sciences)
------------------------------------------------------------------------- */
#include "string.h"
#include "stdlib.h"
#include "compute_cna_atom.h"
#include "atom.h"
#include "update.h"
#include "force.h"
#include "pair.h"
#include "modify.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "neigh_request.h"
#include "comm.h"
#include "memory.h"
#include "error.h"
#include "math.h"
using namespace LAMMPS_NS;
#define MAXNEAR 16
#define MAXCOMMON 8
enum{UNKNOWN,FCC,HCP,BCC,ICOS,OTHER};
enum{NCOMMON,NBOND,MAXBOND,MINBOND};
/* ---------------------------------------------------------------------- */
ComputeCNAAtom::ComputeCNAAtom(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg != 4) error->all(FLERR,"Illegal compute cna/atom command");
peratom_flag = 1;
size_peratom_cols = 0;
double cutoff = atof(arg[3]);
if (cutoff < 0.0) error->all(FLERR,"Illegal compute cna/atom command");
cutsq = cutoff*cutoff;
nmax = 0;
nearest = NULL;
nnearest = NULL;
pattern = NULL;
}
/* ---------------------------------------------------------------------- */
ComputeCNAAtom::~ComputeCNAAtom()
{
memory->destroy(nearest);
memory->destroy(nnearest);
memory->destroy(pattern);
}
/* ---------------------------------------------------------------------- */
void ComputeCNAAtom::init()
{
if (force->pair == NULL)
error->all(FLERR,"Compute cna/atom requires a pair style be defined");
if (sqrt(cutsq) > force->pair->cutforce)
error->all(FLERR,"Compute cna/atom cutoff is longer than pairwise cutoff");
// cannot use neighbor->cutneighmax b/c neighbor has not yet been init
if (2.0*sqrt(cutsq) > force->pair->cutforce + neighbor->skin &&
comm->me == 0)
error->warning(FLERR,"Compute cna/atom cutoff may be too large to find "
"ghost atom neighbors");
int count = 0;
for (int i = 0; i < modify->ncompute; i++)
if (strcmp(modify->compute[i]->style,"cna/atom") == 0) count++;
if (count > 1 && comm->me == 0)
error->warning(FLERR,"More than one compute cna/atom defined");
// need an occasional full neighbor list
int irequest = neighbor->request((void *) this);
neighbor->requests[irequest]->pair = 0;
neighbor->requests[irequest]->compute = 1;
neighbor->requests[irequest]->half = 0;
neighbor->requests[irequest]->full = 1;
neighbor->requests[irequest]->occasional = 1;
}
/* ---------------------------------------------------------------------- */
void ComputeCNAAtom::init_list(int id, NeighList *ptr)
{
list = ptr;
}
/* ---------------------------------------------------------------------- */
void ComputeCNAAtom::compute_peratom()
{
int i,j,k,ii,jj,kk,m,n,inum,jnum,inear,jnear;
int firstflag,ncommon,nbonds,maxbonds,minbonds;
int nfcc,nhcp,nbcc4,nbcc6,nico,cj,ck,cl,cm;
int *ilist,*jlist,*numneigh,**firstneigh;
int cna[MAXNEAR][4],onenearest[MAXNEAR];
int common[MAXCOMMON],bonds[MAXCOMMON];
double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
invoked_peratom = update->ntimestep;
// grow arrays if necessary
if (atom->nlocal > nmax) {
memory->destroy(nearest);
memory->destroy(nnearest);
memory->destroy(pattern);
nmax = atom->nmax;
memory->create(nearest,nmax,MAXNEAR,"cna:nearest");
memory->create(nnearest,nmax,"cna:nnearest");
memory->create(pattern,nmax,"cna:cna_pattern");
vector_atom = pattern;
}
// invoke full neighbor list (will copy or build if necessary)
neighbor->build_one(list->index);
inum = list->inum;
ilist = list->ilist;
numneigh = list->numneigh;
firstneigh = list->firstneigh;
// find the neigbours of each atom within cutoff using full neighbor list
// nearest[] = atom indices of nearest neighbors, up to MAXNEAR
// do this for all atoms, not just compute group
// since CNA calculation requires neighbors of neighbors
double **x = atom->x;
int *mask = atom->mask;
int nlocal = atom->nlocal;
int nerror = 0;
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];
n = 0;
for (jj = 0; jj < jnum; jj++) {
j = jlist[jj];
j &= NEIGHMASK;
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
delz = ztmp - x[j][2];
rsq = delx*delx + dely*dely + delz*delz;
if (rsq < cutsq) {
if (n < MAXNEAR) nearest[i][n++] = j;
else {
nerror++;
break;
}
}
}
nnearest[i] = n;
}
// warning message
int nerrorall;
MPI_Allreduce(&nerror,&nerrorall,1,MPI_INT,MPI_SUM,world);
if (nerrorall && comm->me == 0) {
char str[128];
sprintf(str,"Too many neighbors in CNA for %d atoms",nerrorall);
error->warning(FLERR,str,0);
}
// compute CNA for each atom in group
// only performed if # of nearest neighbors = 12 or 14 (fcc,hcp)
nerror = 0;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
if (!(mask[i] & groupbit)) {
pattern[i] = UNKNOWN;
continue;
}
if (nnearest[i] != 12 && nnearest[i] != 14) {
pattern[i] = OTHER;
continue;
}
// loop over near neighbors of I to build cna data structure
// cna[k][NCOMMON] = # of common neighbors of I with each of its neighs
// cna[k][NBONDS] = # of bonds between those common neighbors
// cna[k][MAXBOND] = max # of bonds of any common neighbor
// cna[k][MINBOND] = min # of bonds of any common neighbor
for (m = 0; m < nnearest[i]; m++) {
j = nearest[i][m];
// common = list of neighbors common to atom I and atom J
// if J is an owned atom, use its near neighbor list to find them
// if J is a ghost atom, use full neighbor list of I to find them
// in latter case, must exclude J from I's neighbor list
if (j < nlocal) {
firstflag = 1;
ncommon = 0;
for (inear = 0; inear < nnearest[i]; inear++)
for (jnear = 0; jnear < nnearest[j]; jnear++)
if (nearest[i][inear] == nearest[j][jnear]) {
if (ncommon < MAXCOMMON) common[ncommon++] = nearest[i][inear];
else if (firstflag) {
nerror++;
firstflag = 0;
}
}
} else {
xtmp = x[j][0];
ytmp = x[j][1];
ztmp = x[j][2];
jlist = firstneigh[i];
jnum = numneigh[i];
n = 0;
for (kk = 0; kk < jnum; kk++) {
k = jlist[kk];
k &= NEIGHMASK;
if (k == j) continue;
delx = xtmp - x[k][0];
dely = ytmp - x[k][1];
delz = ztmp - x[k][2];
rsq = delx*delx + dely*dely + delz*delz;
if (rsq < cutsq) {
if (n < MAXNEAR) onenearest[n++] = k;
else break;
}
}
firstflag = 1;
ncommon = 0;
for (inear = 0; inear < nnearest[i]; inear++)
for (jnear = 0; jnear < n; jnear++)
if (nearest[i][inear] == onenearest[jnear]) {
if (ncommon < MAXCOMMON) common[ncommon++] = nearest[i][inear];
else if (firstflag) {
nerror++;
firstflag = 0;
}
}
}
cna[m][NCOMMON] = ncommon;
// calculate total # of bonds between common neighbor atoms
// also max and min # of common atoms any common atom is bonded to
// bond = pair of atoms within cutoff
for (n = 0; n < ncommon; n++) bonds[n] = 0;
nbonds = 0;
for (jj = 0; jj < ncommon; jj++) {
j = common[jj];
xtmp = x[j][0];
ytmp = x[j][1];
ztmp = x[j][2];
for (kk = jj+1; kk < ncommon; kk++) {
k = common[kk];
delx = xtmp - x[k][0];
dely = ytmp - x[k][1];
delz = ztmp - x[k][2];
rsq = delx*delx + dely*dely + delz*delz;
if (rsq < cutsq) {
nbonds++;
bonds[jj]++;
bonds[kk]++;
}
}
}
cna[m][NBOND] = nbonds;
maxbonds = 0;
minbonds = MAXCOMMON;
for (n = 0; n < ncommon; n++) {
maxbonds = MAX(bonds[n],maxbonds);
minbonds = MIN(bonds[n],minbonds);
}
cna[m][MAXBOND] = maxbonds;
cna[m][MINBOND] = minbonds;
}
// detect CNA pattern of the atom
nfcc = nhcp = nbcc4 = nbcc6 = nico = 0;
pattern[i] = OTHER;
if (nnearest[i] == 12) {
for (inear = 0; inear < 12; inear++) {
cj = cna[inear][NCOMMON];
ck = cna[inear][NBOND];
cl = cna[inear][MAXBOND];
cm = cna[inear][MINBOND];
if (cj == 4 && ck == 2 && cl == 1 && cm == 1) nfcc++;
else if (cj == 4 && ck == 2 && cl == 2 && cm == 0) nhcp++;
else if (cj == 5 && ck == 5 && cl == 2 && cm == 2) nico++;
}
if (nfcc == 12) pattern[i] = FCC;
else if (nfcc == 6 && nhcp == 6) pattern[i] = HCP;
else if (nico == 12) pattern[i] = ICOS;
} else if (nnearest[i] == 14) {
for (inear = 0; inear < 14; inear++) {
cj = cna[inear][NCOMMON];
ck = cna[inear][NBOND];
cl = cna[inear][MAXBOND];
cm = cna[inear][MINBOND];
if (cj == 4 && ck == 4 && cl == 2 && cm == 2) nbcc4++;
else if (cj == 6 && ck == 6 && cl == 2 && cm == 2) nbcc6++;
}
if (nbcc4 == 6 && nbcc6 == 8) pattern[i] = BCC;
}
}
// warning message
MPI_Allreduce(&nerror,&nerrorall,1,MPI_INT,MPI_SUM,world);
if (nerrorall && comm->me == 0) {
char str[128];
sprintf(str,"Too many common neighbors in CNA %d times",nerrorall);
error->warning(FLERR,str);
}
}
/* ----------------------------------------------------------------------
memory usage of local atom-based array
------------------------------------------------------------------------- */
double ComputeCNAAtom::memory_usage()
{
double bytes = nmax * sizeof(int);
bytes += nmax * MAXNEAR * sizeof(int);
bytes += nmax * sizeof(double);
return bytes;
}
|