convert more styles to use new neighbor list request API, apply clang-format
This commit is contained in:
@ -1,4 +1,3 @@
|
||||
// clang-format off
|
||||
/* ----------------------------------------------------------------------
|
||||
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
||||
https://www.lammps.org/, Sandia National Laboratories
|
||||
@ -28,23 +27,21 @@
|
||||
#include "memory.h"
|
||||
#include "modify.h"
|
||||
#include "neigh_list.h"
|
||||
#include "neigh_request.h"
|
||||
#include "neighbor.h"
|
||||
#include "pair.h"
|
||||
#include "update.h"
|
||||
|
||||
#include <cstring>
|
||||
#include <cmath>
|
||||
#include <cstring>
|
||||
|
||||
using namespace LAMMPS_NS;
|
||||
using namespace MathConst;
|
||||
using namespace MathSpecial;
|
||||
|
||||
|
||||
#ifdef DBL_EPSILON
|
||||
#define MY_EPSILON (10.0*DBL_EPSILON)
|
||||
#define MY_EPSILON (10.0 * DBL_EPSILON)
|
||||
#else
|
||||
#define MY_EPSILON (10.0*2.220446049250313e-16)
|
||||
#define MY_EPSILON (10.0 * 2.220446049250313e-16)
|
||||
#endif
|
||||
|
||||
#define QEPSILON 1.0e-6
|
||||
@ -52,11 +49,10 @@ using namespace MathSpecial;
|
||||
/* ---------------------------------------------------------------------- */
|
||||
|
||||
ComputeOrientOrderAtom::ComputeOrientOrderAtom(LAMMPS *lmp, int narg, char **arg) :
|
||||
Compute(lmp, narg, arg),
|
||||
qlist(nullptr), distsq(nullptr), nearest(nullptr), rlist(nullptr),
|
||||
qnarray(nullptr), qnm_r(nullptr), qnm_i(nullptr), cglist(nullptr)
|
||||
Compute(lmp, narg, arg), qlist(nullptr), distsq(nullptr), nearest(nullptr), rlist(nullptr),
|
||||
qnarray(nullptr), qnm_r(nullptr), qnm_i(nullptr), cglist(nullptr)
|
||||
{
|
||||
if (narg < 3 ) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
if (narg < 3) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
|
||||
// set default values for optional args
|
||||
|
||||
@ -70,7 +66,7 @@ ComputeOrientOrderAtom::ComputeOrientOrderAtom(LAMMPS *lmp, int narg, char **arg
|
||||
// specify which orders to request
|
||||
|
||||
nqlist = 5;
|
||||
memory->create(qlist,nqlist,"orientorder/atom:qlist");
|
||||
memory->create(qlist, nqlist, "orientorder/atom:qlist");
|
||||
qlist[0] = 4;
|
||||
qlist[1] = 6;
|
||||
qlist[2] = 8;
|
||||
@ -82,69 +78,69 @@ ComputeOrientOrderAtom::ComputeOrientOrderAtom(LAMMPS *lmp, int narg, char **arg
|
||||
|
||||
int iarg = 3;
|
||||
while (iarg < narg) {
|
||||
if (strcmp(arg[iarg],"nnn") == 0) {
|
||||
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
if (strcmp(arg[iarg+1],"NULL") == 0) {
|
||||
if (strcmp(arg[iarg], "nnn") == 0) {
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
if (strcmp(arg[iarg + 1], "NULL") == 0) {
|
||||
nnn = 0;
|
||||
} else {
|
||||
nnn = utils::numeric(FLERR,arg[iarg+1],false,lmp);
|
||||
if (nnn <= 0) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
nnn = utils::numeric(FLERR, arg[iarg + 1], false, lmp);
|
||||
if (nnn <= 0) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
}
|
||||
iarg += 2;
|
||||
} else if (strcmp(arg[iarg],"degrees") == 0) {
|
||||
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
nqlist = utils::numeric(FLERR,arg[iarg+1],false,lmp);
|
||||
if (nqlist <= 0) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
} else if (strcmp(arg[iarg], "degrees") == 0) {
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
nqlist = utils::numeric(FLERR, arg[iarg + 1], false, lmp);
|
||||
if (nqlist <= 0) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
memory->destroy(qlist);
|
||||
memory->create(qlist,nqlist,"orientorder/atom:qlist");
|
||||
memory->create(qlist, nqlist, "orientorder/atom:qlist");
|
||||
iarg += 2;
|
||||
if (iarg+nqlist > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
if (iarg + nqlist > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
qmax = 0;
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
qlist[il] = utils::numeric(FLERR,arg[iarg+il],false,lmp);
|
||||
if (qlist[il] < 0) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
qlist[il] = utils::numeric(FLERR, arg[iarg + il], false, lmp);
|
||||
if (qlist[il] < 0) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
if (qlist[il] > qmax) qmax = qlist[il];
|
||||
}
|
||||
iarg += nqlist;
|
||||
} else if (strcmp(arg[iarg],"wl") == 0) {
|
||||
if (iarg+2 > narg)
|
||||
error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
wlflag = utils::logical(FLERR,arg[iarg+1],false,lmp);
|
||||
} else if (strcmp(arg[iarg], "wl") == 0) {
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
wlflag = utils::logical(FLERR, arg[iarg + 1], false, lmp);
|
||||
iarg += 2;
|
||||
} else if (strcmp(arg[iarg],"wl/hat") == 0) {
|
||||
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
wlhatflag = utils::logical(FLERR,arg[iarg+1],false,lmp);
|
||||
} else if (strcmp(arg[iarg], "wl/hat") == 0) {
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
wlhatflag = utils::logical(FLERR, arg[iarg + 1], false, lmp);
|
||||
iarg += 2;
|
||||
} else if (strcmp(arg[iarg],"components") == 0) {
|
||||
} else if (strcmp(arg[iarg], "components") == 0) {
|
||||
qlcompflag = 1;
|
||||
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
qlcomp = utils::numeric(FLERR,arg[iarg+1],false,lmp);
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
qlcomp = utils::numeric(FLERR, arg[iarg + 1], false, lmp);
|
||||
iqlcomp = -1;
|
||||
for (int il = 0; il < nqlist; il++)
|
||||
if (qlcomp == qlist[il]) {
|
||||
iqlcomp = il;
|
||||
break;
|
||||
}
|
||||
if (iqlcomp == -1) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
if (iqlcomp == -1) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
iarg += 2;
|
||||
} else if (strcmp(arg[iarg],"cutoff") == 0) {
|
||||
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
double cutoff = utils::numeric(FLERR,arg[iarg+1],false,lmp);
|
||||
if (cutoff <= 0.0) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
cutsq = cutoff*cutoff;
|
||||
} else if (strcmp(arg[iarg], "cutoff") == 0) {
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
double cutoff = utils::numeric(FLERR, arg[iarg + 1], false, lmp);
|
||||
if (cutoff <= 0.0) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
cutsq = cutoff * cutoff;
|
||||
iarg += 2;
|
||||
} else if (strcmp(arg[iarg],"chunksize") == 0) {
|
||||
if (iarg+2 > narg) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
chunksize = utils::numeric(FLERR,arg[iarg+1],false,lmp);
|
||||
if (chunksize <= 0) error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
} else if (strcmp(arg[iarg], "chunksize") == 0) {
|
||||
if (iarg + 2 > narg) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
chunksize = utils::numeric(FLERR, arg[iarg + 1], false, lmp);
|
||||
if (chunksize <= 0) error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
iarg += 2;
|
||||
} else error->all(FLERR,"Illegal compute orientorder/atom command");
|
||||
} else
|
||||
error->all(FLERR, "Illegal compute orientorder/atom command");
|
||||
}
|
||||
|
||||
ncol = nqlist;
|
||||
if (wlflag) ncol += nqlist;
|
||||
if (wlhatflag) ncol += nqlist;
|
||||
if (qlcompflag) ncol += 2*(2*qlcomp+1);
|
||||
if (qlcompflag) ncol += 2 * (2 * qlcomp + 1);
|
||||
|
||||
peratom_flag = 1;
|
||||
size_peratom_cols = ncol;
|
||||
@ -174,32 +170,23 @@ ComputeOrientOrderAtom::~ComputeOrientOrderAtom()
|
||||
void ComputeOrientOrderAtom::init()
|
||||
{
|
||||
if (force->pair == nullptr)
|
||||
error->all(FLERR,"Compute orientorder/atom requires a "
|
||||
"pair style be defined");
|
||||
if (cutsq == 0.0) cutsq = force->pair->cutforce * force->pair->cutforce;
|
||||
error->all(FLERR, "Compute orientorder/atom requires a pair style be defined");
|
||||
if (cutsq == 0.0)
|
||||
cutsq = force->pair->cutforce * force->pair->cutforce;
|
||||
else if (sqrt(cutsq) > force->pair->cutforce)
|
||||
error->all(FLERR,"Compute orientorder/atom cutoff is "
|
||||
"longer than pairwise cutoff");
|
||||
error->all(FLERR, "Compute orientorder/atom cutoff is longer than pairwise cutoff");
|
||||
|
||||
memory->destroy(qnm_r);
|
||||
memory->destroy(qnm_i);
|
||||
memory->create(qnm_r,nqlist,2*qmax+1,"orientorder/atom:qnm_r");
|
||||
memory->create(qnm_i,nqlist,2*qmax+1,"orientorder/atom:qnm_i");
|
||||
memory->create(qnm_r, nqlist, 2 * qmax + 1, "orientorder/atom:qnm_r");
|
||||
memory->create(qnm_i, nqlist, 2 * qmax + 1, "orientorder/atom:qnm_i");
|
||||
|
||||
// need an occasional full neighbor list
|
||||
|
||||
int irequest = neighbor->request(this,instance_me);
|
||||
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;
|
||||
neighbor->add_request(this, NeighConst::REQ_FULL | NeighConst::REQ_OCCASIONAL);
|
||||
|
||||
int count = 0;
|
||||
for (int i = 0; i < modify->ncompute; i++)
|
||||
if (strcmp(modify->compute[i]->style,"orientorder/atom") == 0) count++;
|
||||
if (count > 1 && comm->me == 0)
|
||||
error->warning(FLERR,"More than one compute orientorder/atom");
|
||||
if ((modify->get_compute_by_style("orientorder/atom").size() > 1) && (comm->me == 0))
|
||||
error->warning(FLERR, "More than one instance of compute orientorder/atom");
|
||||
|
||||
if (wlflag || wlhatflag) init_clebsch_gordan();
|
||||
}
|
||||
@ -215,9 +202,9 @@ void ComputeOrientOrderAtom::init_list(int /*id*/, NeighList *ptr)
|
||||
|
||||
void ComputeOrientOrderAtom::compute_peratom()
|
||||
{
|
||||
int i,j,ii,jj,inum,jnum;
|
||||
double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
|
||||
int *ilist,*jlist,*numneigh,**firstneigh;
|
||||
int i, j, ii, jj, inum, jnum;
|
||||
double xtmp, ytmp, ztmp, delx, dely, delz, rsq;
|
||||
int *ilist, *jlist, *numneigh, **firstneigh;
|
||||
|
||||
invoked_peratom = update->ntimestep;
|
||||
|
||||
@ -226,7 +213,7 @@ void ComputeOrientOrderAtom::compute_peratom()
|
||||
if (atom->nmax > nmax) {
|
||||
memory->destroy(qnarray);
|
||||
nmax = atom->nmax;
|
||||
memory->create(qnarray,nmax,ncol,"orientorder/atom:qnarray");
|
||||
memory->create(qnarray, nmax, ncol, "orientorder/atom:qnarray");
|
||||
array_atom = qnarray;
|
||||
}
|
||||
|
||||
@ -244,11 +231,11 @@ void ComputeOrientOrderAtom::compute_peratom()
|
||||
|
||||
double **x = atom->x;
|
||||
int *mask = atom->mask;
|
||||
memset(&qnarray[0][0],0,sizeof(double)*nmax*ncol);
|
||||
memset(&qnarray[0][0], 0, sizeof(double) * nmax * ncol);
|
||||
|
||||
for (ii = 0; ii < inum; ii++) {
|
||||
i = ilist[ii];
|
||||
double* qn = qnarray[i];
|
||||
double *qn = qnarray[i];
|
||||
if (mask[i] & groupbit) {
|
||||
xtmp = x[i][0];
|
||||
ytmp = x[i][1];
|
||||
@ -263,9 +250,9 @@ void ComputeOrientOrderAtom::compute_peratom()
|
||||
memory->destroy(rlist);
|
||||
memory->destroy(nearest);
|
||||
maxneigh = jnum;
|
||||
memory->create(distsq,maxneigh,"orientorder/atom:distsq");
|
||||
memory->create(rlist,maxneigh,3,"orientorder/atom:rlist");
|
||||
memory->create(nearest,maxneigh,"orientorder/atom:nearest");
|
||||
memory->create(distsq, maxneigh, "orientorder/atom:distsq");
|
||||
memory->create(rlist, maxneigh, 3, "orientorder/atom:rlist");
|
||||
memory->create(nearest, maxneigh, "orientorder/atom:nearest");
|
||||
}
|
||||
|
||||
// loop over list of all neighbors within force cutoff
|
||||
@ -281,7 +268,7 @@ void ComputeOrientOrderAtom::compute_peratom()
|
||||
delx = xtmp - x[j][0];
|
||||
dely = ytmp - x[j][1];
|
||||
delz = ztmp - x[j][2];
|
||||
rsq = delx*delx + dely*dely + delz*delz;
|
||||
rsq = delx * delx + dely * dely + delz * delz;
|
||||
if (rsq < cutsq) {
|
||||
distsq[ncount] = rsq;
|
||||
rlist[ncount][0] = delx;
|
||||
@ -294,15 +281,14 @@ void ComputeOrientOrderAtom::compute_peratom()
|
||||
// if not nnn neighbors, order parameter = 0;
|
||||
|
||||
if ((ncount == 0) || (ncount < nnn)) {
|
||||
for (jj = 0; jj < ncol; jj++)
|
||||
qn[jj] = 0.0;
|
||||
for (jj = 0; jj < ncol; jj++) qn[jj] = 0.0;
|
||||
continue;
|
||||
}
|
||||
|
||||
// if nnn > 0, use only nearest nnn neighbors
|
||||
|
||||
if (nnn > 0) {
|
||||
select3(nnn,ncount,distsq,nearest,rlist);
|
||||
select3(nnn, ncount, distsq, nearest, rlist);
|
||||
ncount = nnn;
|
||||
}
|
||||
|
||||
@ -317,9 +303,9 @@ void ComputeOrientOrderAtom::compute_peratom()
|
||||
|
||||
double ComputeOrientOrderAtom::memory_usage()
|
||||
{
|
||||
double bytes = (double)ncol*nmax * sizeof(double);
|
||||
bytes += (double)(qmax*(2*qmax+1)+maxneigh*4) * sizeof(double);
|
||||
bytes += (double)(nqlist+maxneigh) * sizeof(int);
|
||||
double bytes = (double) ncol * nmax * sizeof(double);
|
||||
bytes += (double) (qmax * (2 * qmax + 1) + maxneigh * 4) * sizeof(double);
|
||||
bytes += (double) (nqlist + maxneigh) * sizeof(int);
|
||||
return bytes;
|
||||
}
|
||||
|
||||
@ -331,26 +317,39 @@ double ComputeOrientOrderAtom::memory_usage()
|
||||
|
||||
// Use no-op do while to create single statement
|
||||
|
||||
#define SWAP(a,b) do { \
|
||||
tmp = a; a = b; b = tmp; \
|
||||
#define SWAP(a, b) \
|
||||
do { \
|
||||
tmp = a; \
|
||||
a = b; \
|
||||
b = tmp; \
|
||||
} while (0)
|
||||
|
||||
#define ISWAP(a,b) do { \
|
||||
itmp = a; a = b; b = itmp; \
|
||||
#define ISWAP(a, b) \
|
||||
do { \
|
||||
itmp = a; \
|
||||
a = b; \
|
||||
b = itmp; \
|
||||
} while (0)
|
||||
|
||||
#define SWAP3(a,b) do { \
|
||||
tmp = a[0]; a[0] = b[0]; b[0] = tmp; \
|
||||
tmp = a[1]; a[1] = b[1]; b[1] = tmp; \
|
||||
tmp = a[2]; a[2] = b[2]; b[2] = tmp; \
|
||||
#define SWAP3(a, b) \
|
||||
do { \
|
||||
tmp = a[0]; \
|
||||
a[0] = b[0]; \
|
||||
b[0] = tmp; \
|
||||
tmp = a[1]; \
|
||||
a[1] = b[1]; \
|
||||
b[1] = tmp; \
|
||||
tmp = a[2]; \
|
||||
a[2] = b[2]; \
|
||||
b[2] = tmp; \
|
||||
} while (0)
|
||||
|
||||
/* ---------------------------------------------------------------------- */
|
||||
|
||||
void ComputeOrientOrderAtom::select3(int k, int n, double *arr, int *iarr, double **arr3)
|
||||
{
|
||||
int i,ir,j,l,mid,ia,itmp;
|
||||
double a,tmp,a3[3];
|
||||
int i, ir, j, l, mid, ia, itmp;
|
||||
double a, tmp, a3[3];
|
||||
|
||||
arr--;
|
||||
iarr--;
|
||||
@ -358,59 +357,61 @@ void ComputeOrientOrderAtom::select3(int k, int n, double *arr, int *iarr, doubl
|
||||
l = 1;
|
||||
ir = n;
|
||||
for (;;) {
|
||||
if (ir <= l+1) {
|
||||
if (ir == l+1 && arr[ir] < arr[l]) {
|
||||
SWAP(arr[l],arr[ir]);
|
||||
ISWAP(iarr[l],iarr[ir]);
|
||||
SWAP3(arr3[l],arr3[ir]);
|
||||
if (ir <= l + 1) {
|
||||
if (ir == l + 1 && arr[ir] < arr[l]) {
|
||||
SWAP(arr[l], arr[ir]);
|
||||
ISWAP(iarr[l], iarr[ir]);
|
||||
SWAP3(arr3[l], arr3[ir]);
|
||||
}
|
||||
return;
|
||||
} else {
|
||||
mid=(l+ir) >> 1;
|
||||
SWAP(arr[mid],arr[l+1]);
|
||||
ISWAP(iarr[mid],iarr[l+1]);
|
||||
SWAP3(arr3[mid],arr3[l+1]);
|
||||
mid = (l + ir) >> 1;
|
||||
SWAP(arr[mid], arr[l + 1]);
|
||||
ISWAP(iarr[mid], iarr[l + 1]);
|
||||
SWAP3(arr3[mid], arr3[l + 1]);
|
||||
if (arr[l] > arr[ir]) {
|
||||
SWAP(arr[l],arr[ir]);
|
||||
ISWAP(iarr[l],iarr[ir]);
|
||||
SWAP3(arr3[l],arr3[ir]);
|
||||
SWAP(arr[l], arr[ir]);
|
||||
ISWAP(iarr[l], iarr[ir]);
|
||||
SWAP3(arr3[l], arr3[ir]);
|
||||
}
|
||||
if (arr[l+1] > arr[ir]) {
|
||||
SWAP(arr[l+1],arr[ir]);
|
||||
ISWAP(iarr[l+1],iarr[ir]);
|
||||
SWAP3(arr3[l+1],arr3[ir]);
|
||||
if (arr[l + 1] > arr[ir]) {
|
||||
SWAP(arr[l + 1], arr[ir]);
|
||||
ISWAP(iarr[l + 1], iarr[ir]);
|
||||
SWAP3(arr3[l + 1], arr3[ir]);
|
||||
}
|
||||
if (arr[l] > arr[l+1]) {
|
||||
SWAP(arr[l],arr[l+1]);
|
||||
ISWAP(iarr[l],iarr[l+1]);
|
||||
SWAP3(arr3[l],arr3[l+1]);
|
||||
if (arr[l] > arr[l + 1]) {
|
||||
SWAP(arr[l], arr[l + 1]);
|
||||
ISWAP(iarr[l], iarr[l + 1]);
|
||||
SWAP3(arr3[l], arr3[l + 1]);
|
||||
}
|
||||
i = l+1;
|
||||
i = l + 1;
|
||||
j = ir;
|
||||
a = arr[l+1];
|
||||
ia = iarr[l+1];
|
||||
a3[0] = arr3[l+1][0];
|
||||
a3[1] = arr3[l+1][1];
|
||||
a3[2] = arr3[l+1][2];
|
||||
a = arr[l + 1];
|
||||
ia = iarr[l + 1];
|
||||
a3[0] = arr3[l + 1][0];
|
||||
a3[1] = arr3[l + 1][1];
|
||||
a3[2] = arr3[l + 1][2];
|
||||
for (;;) {
|
||||
do i++; while (arr[i] < a);
|
||||
do j--; while (arr[j] > a);
|
||||
do i++;
|
||||
while (arr[i] < a);
|
||||
do j--;
|
||||
while (arr[j] > a);
|
||||
if (j < i) break;
|
||||
SWAP(arr[i],arr[j]);
|
||||
ISWAP(iarr[i],iarr[j]);
|
||||
SWAP3(arr3[i],arr3[j]);
|
||||
SWAP(arr[i], arr[j]);
|
||||
ISWAP(iarr[i], iarr[j]);
|
||||
SWAP3(arr3[i], arr3[j]);
|
||||
}
|
||||
arr[l+1] = arr[j];
|
||||
arr[l + 1] = arr[j];
|
||||
arr[j] = a;
|
||||
iarr[l+1] = iarr[j];
|
||||
iarr[l + 1] = iarr[j];
|
||||
iarr[j] = ia;
|
||||
arr3[l+1][0] = arr3[j][0];
|
||||
arr3[l+1][1] = arr3[j][1];
|
||||
arr3[l+1][2] = arr3[j][2];
|
||||
arr3[l + 1][0] = arr3[j][0];
|
||||
arr3[l + 1][1] = arr3[j][1];
|
||||
arr3[l + 1][2] = arr3[j][2];
|
||||
arr3[j][0] = a3[0];
|
||||
arr3[j][1] = a3[1];
|
||||
arr3[j][2] = a3[2];
|
||||
if (j >= k) ir = j-1;
|
||||
if (j >= k) ir = j - 1;
|
||||
if (j <= k) l = i;
|
||||
}
|
||||
}
|
||||
@ -420,34 +421,32 @@ void ComputeOrientOrderAtom::select3(int k, int n, double *arr, int *iarr, doubl
|
||||
calculate the bond orientational order parameters
|
||||
------------------------------------------------------------------------- */
|
||||
|
||||
void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
int ncount, double qn[],
|
||||
int qlist[], int nqlist) {
|
||||
void ComputeOrientOrderAtom::calc_boop(double **rlist, int ncount, double qn[], int qlist[],
|
||||
int nqlist)
|
||||
{
|
||||
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
for (int m = 0; m < 2*l+1; m++) {
|
||||
for (int m = 0; m < 2 * l + 1; m++) {
|
||||
qnm_r[il][m] = 0.0;
|
||||
qnm_i[il][m] = 0.0;
|
||||
}
|
||||
}
|
||||
|
||||
for (int ineigh = 0; ineigh < ncount; ineigh++) {
|
||||
const double * const r = rlist[ineigh];
|
||||
double rmag = dist(r);
|
||||
if (rmag <= MY_EPSILON) {
|
||||
return;
|
||||
}
|
||||
const double *const r = rlist[ineigh];
|
||||
double rmag = sqrt(r[0] * r[0] + r[1] * r[1] + r[2] * r[2]);
|
||||
if (rmag <= MY_EPSILON) { return; }
|
||||
|
||||
double costheta = r[2] / rmag;
|
||||
double expphi_r = r[0];
|
||||
double expphi_i = r[1];
|
||||
double rxymag = sqrt(expphi_r*expphi_r+expphi_i*expphi_i);
|
||||
double rxymag = sqrt(expphi_r * expphi_r + expphi_i * expphi_i);
|
||||
if (rxymag <= MY_EPSILON) {
|
||||
expphi_r = 1.0;
|
||||
expphi_i = 0.0;
|
||||
} else {
|
||||
double rxymaginv = 1.0/rxymag;
|
||||
double rxymaginv = 1.0 / rxymag;
|
||||
expphi_r *= rxymaginv;
|
||||
expphi_i *= rxymaginv;
|
||||
}
|
||||
@ -467,21 +466,20 @@ void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
double prefactor = polar_prefactor(l, m, costheta);
|
||||
double ylm_r = prefactor * expphim_r;
|
||||
double ylm_i = prefactor * expphim_i;
|
||||
qnm_r[il][m+l] += ylm_r;
|
||||
qnm_i[il][m+l] += ylm_i;
|
||||
qnm_r[il][m + l] += ylm_r;
|
||||
qnm_i[il][m + l] += ylm_i;
|
||||
if (m & 1) {
|
||||
qnm_r[il][-m+l] -= ylm_r;
|
||||
qnm_i[il][-m+l] += ylm_i;
|
||||
qnm_r[il][-m + l] -= ylm_r;
|
||||
qnm_i[il][-m + l] += ylm_i;
|
||||
} else {
|
||||
qnm_r[il][-m+l] += ylm_r;
|
||||
qnm_i[il][-m+l] -= ylm_i;
|
||||
qnm_r[il][-m + l] += ylm_r;
|
||||
qnm_i[il][-m + l] -= ylm_i;
|
||||
}
|
||||
double tmp_r = expphim_r*expphi_r - expphim_i*expphi_i;
|
||||
double tmp_i = expphim_r*expphi_i + expphim_i*expphi_r;
|
||||
double tmp_r = expphim_r * expphi_r - expphim_i * expphi_i;
|
||||
double tmp_i = expphim_r * expphi_i + expphim_i * expphi_r;
|
||||
expphim_r = tmp_r;
|
||||
expphim_i = tmp_i;
|
||||
}
|
||||
|
||||
}
|
||||
}
|
||||
|
||||
@ -490,7 +488,7 @@ void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
double facn = 1.0 / ncount;
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
for (int m = 0; m < 2*l+1; m++) {
|
||||
for (int m = 0; m < 2 * l + 1; m++) {
|
||||
qnm_r[il][m] *= facn;
|
||||
qnm_i[il][m] *= facn;
|
||||
}
|
||||
@ -502,10 +500,10 @@ void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
int jj = 0;
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
double qnormfac = sqrt(MY_4PI/(2*l+1));
|
||||
double qnormfac = sqrt(MY_4PI / (2 * l + 1));
|
||||
double qm_sum = 0.0;
|
||||
for (int m = 0; m < 2*l+1; m++)
|
||||
qm_sum += qnm_r[il][m]*qnm_r[il][m] + qnm_i[il][m]*qnm_i[il][m];
|
||||
for (int m = 0; m < 2 * l + 1; m++)
|
||||
qm_sum += qnm_r[il][m] * qnm_r[il][m] + qnm_i[il][m] * qnm_i[il][m];
|
||||
qn[jj++] = qnormfac * sqrt(qm_sum);
|
||||
}
|
||||
|
||||
@ -516,16 +514,16 @@ void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
double wlsum = 0.0;
|
||||
for (int m1 = 0; m1 < 2*l+1; m1++) {
|
||||
for (int m2 = MAX(0,l-m1); m2 < MIN(2*l+1,3*l-m1+1); m2++) {
|
||||
for (int m1 = 0; m1 < 2 * l + 1; m1++) {
|
||||
for (int m2 = MAX(0, l - m1); m2 < MIN(2 * l + 1, 3 * l - m1 + 1); m2++) {
|
||||
int m = m1 + m2 - l;
|
||||
double qm1qm2_r = qnm_r[il][m1]*qnm_r[il][m2] - qnm_i[il][m1]*qnm_i[il][m2];
|
||||
double qm1qm2_i = qnm_r[il][m1]*qnm_i[il][m2] + qnm_i[il][m1]*qnm_r[il][m2];
|
||||
wlsum += (qm1qm2_r*qnm_r[il][m] + qm1qm2_i*qnm_i[il][m])*cglist[idxcg_count];
|
||||
double qm1qm2_r = qnm_r[il][m1] * qnm_r[il][m2] - qnm_i[il][m1] * qnm_i[il][m2];
|
||||
double qm1qm2_i = qnm_r[il][m1] * qnm_i[il][m2] + qnm_i[il][m1] * qnm_r[il][m2];
|
||||
wlsum += (qm1qm2_r * qnm_r[il][m] + qm1qm2_i * qnm_i[il][m]) * cglist[idxcg_count];
|
||||
idxcg_count++;
|
||||
}
|
||||
}
|
||||
qn[jj++] = wlsum/sqrt(2*l+1);
|
||||
qn[jj++] = wlsum / sqrt(2 * l + 1);
|
||||
}
|
||||
}
|
||||
|
||||
@ -536,21 +534,21 @@ void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
double wlsum = 0.0;
|
||||
for (int m1 = 0; m1 < 2*l+1; m1++) {
|
||||
for (int m2 = MAX(0,l-m1); m2 < MIN(2*l+1,3*l-m1+1); m2++) {
|
||||
for (int m1 = 0; m1 < 2 * l + 1; m1++) {
|
||||
for (int m2 = MAX(0, l - m1); m2 < MIN(2 * l + 1, 3 * l - m1 + 1); m2++) {
|
||||
int m = m1 + m2 - l;
|
||||
double qm1qm2_r = qnm_r[il][m1]*qnm_r[il][m2] - qnm_i[il][m1]*qnm_i[il][m2];
|
||||
double qm1qm2_i = qnm_r[il][m1]*qnm_i[il][m2] + qnm_i[il][m1]*qnm_r[il][m2];
|
||||
wlsum += (qm1qm2_r*qnm_r[il][m] + qm1qm2_i*qnm_i[il][m])*cglist[idxcg_count];
|
||||
double qm1qm2_r = qnm_r[il][m1] * qnm_r[il][m2] - qnm_i[il][m1] * qnm_i[il][m2];
|
||||
double qm1qm2_i = qnm_r[il][m1] * qnm_i[il][m2] + qnm_i[il][m1] * qnm_r[il][m2];
|
||||
wlsum += (qm1qm2_r * qnm_r[il][m] + qm1qm2_i * qnm_i[il][m]) * cglist[idxcg_count];
|
||||
idxcg_count++;
|
||||
}
|
||||
}
|
||||
if (qn[il] < QEPSILON)
|
||||
qn[jj++] = 0.0;
|
||||
else {
|
||||
double qnormfac = sqrt(MY_4PI/(2*l+1));
|
||||
double qnfac = qnormfac/qn[il];
|
||||
qn[jj++] = wlsum/sqrt(2*l+1)*(qnfac*qnfac*qnfac);
|
||||
double qnormfac = sqrt(MY_4PI / (2 * l + 1));
|
||||
double qnfac = qnormfac / qn[il];
|
||||
qn[jj++] = wlsum / sqrt(2 * l + 1) * (qnfac * qnfac * qnfac);
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -561,29 +559,19 @@ void ComputeOrientOrderAtom::calc_boop(double **rlist,
|
||||
int il = iqlcomp;
|
||||
int l = qlcomp;
|
||||
if (qn[il] < QEPSILON)
|
||||
for (int m = 0; m < 2*l+1; m++) {
|
||||
for (int m = 0; m < 2 * l + 1; m++) {
|
||||
qn[jj++] = 0.0;
|
||||
qn[jj++] = 0.0;
|
||||
}
|
||||
else {
|
||||
double qnormfac = sqrt(MY_4PI/(2*l+1));
|
||||
double qnfac = qnormfac/qn[il];
|
||||
for (int m = 0; m < 2*l+1; m++) {
|
||||
double qnormfac = sqrt(MY_4PI / (2 * l + 1));
|
||||
double qnfac = qnormfac / qn[il];
|
||||
for (int m = 0; m < 2 * l + 1; m++) {
|
||||
qn[jj++] = qnm_r[il][m] * qnfac;
|
||||
qn[jj++] = qnm_i[il][m] * qnfac;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
/* ----------------------------------------------------------------------
|
||||
calculate scalar distance
|
||||
------------------------------------------------------------------------- */
|
||||
|
||||
double ComputeOrientOrderAtom::dist(const double r[])
|
||||
{
|
||||
return sqrt(r[0]*r[0] + r[1]*r[1] + r[2]*r[2]);
|
||||
}
|
||||
|
||||
/* ----------------------------------------------------------------------
|
||||
@ -596,11 +584,10 @@ double ComputeOrientOrderAtom::polar_prefactor(int l, int m, double costheta)
|
||||
const int mabs = abs(m);
|
||||
|
||||
double prefactor = 1.0;
|
||||
for (int i=l-mabs+1; i < l+mabs+1; ++i)
|
||||
prefactor *= static_cast<double>(i);
|
||||
for (int i = l - mabs + 1; i < l + mabs + 1; ++i) prefactor *= static_cast<double>(i);
|
||||
|
||||
prefactor = sqrt(static_cast<double>(2*l+1)/(MY_4PI*prefactor))
|
||||
* associated_legendre(l,mabs,costheta);
|
||||
prefactor = sqrt(static_cast<double>(2 * l + 1) / (MY_4PI * prefactor)) *
|
||||
associated_legendre(l, mabs, costheta);
|
||||
|
||||
if ((m < 0) && (m % 2)) prefactor = -prefactor;
|
||||
|
||||
@ -619,16 +606,15 @@ double ComputeOrientOrderAtom::associated_legendre(int l, int m, double x)
|
||||
double p(1.0), pm1(0.0), pm2(0.0);
|
||||
|
||||
if (m != 0) {
|
||||
const double msqx = -sqrt(1.0-x*x);
|
||||
for (int i=1; i < m+1; ++i)
|
||||
p *= static_cast<double>(2*i-1) * msqx;
|
||||
const double msqx = -sqrt(1.0 - x * x);
|
||||
for (int i = 1; i < m + 1; ++i) p *= static_cast<double>(2 * i - 1) * msqx;
|
||||
}
|
||||
|
||||
for (int i=m+1; i < l+1; ++i) {
|
||||
for (int i = m + 1; i < l + 1; ++i) {
|
||||
pm2 = pm1;
|
||||
pm1 = p;
|
||||
p = (static_cast<double>(2*i-1)*x*pm1
|
||||
- static_cast<double>(i+m-1)*pm2) / static_cast<double>(i-m);
|
||||
p = (static_cast<double>(2 * i - 1) * x * pm1 - static_cast<double>(i + m - 1) * pm2) /
|
||||
static_cast<double>(i - m);
|
||||
}
|
||||
|
||||
return p;
|
||||
@ -642,16 +628,15 @@ double ComputeOrientOrderAtom::associated_legendre(int l, int m, double x)
|
||||
|
||||
void ComputeOrientOrderAtom::init_clebsch_gordan()
|
||||
{
|
||||
double sum,dcg,sfaccg, sfac1, sfac2;
|
||||
double sum, dcg, sfaccg, sfac1, sfac2;
|
||||
int m, aa2, bb2, cc2;
|
||||
int ifac, idxcg_count;
|
||||
|
||||
idxcg_count = 0;
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
for (int m1 = 0; m1 < 2*l+1; m1++)
|
||||
for (int m2 = MAX(0,l-m1); m2 < MIN(2*l+1,3*l-m1+1); m2++)
|
||||
idxcg_count++;
|
||||
for (int m1 = 0; m1 < 2 * l + 1; m1++)
|
||||
for (int m2 = MAX(0, l - m1); m2 < MIN(2 * l + 1, 3 * l - m1 + 1); m2++) idxcg_count++;
|
||||
}
|
||||
idxcg_max = idxcg_count;
|
||||
memory->destroy(cglist);
|
||||
@ -660,41 +645,44 @@ void ComputeOrientOrderAtom::init_clebsch_gordan()
|
||||
idxcg_count = 0;
|
||||
for (int il = 0; il < nqlist; il++) {
|
||||
int l = qlist[il];
|
||||
for (int m1 = 0; m1 < 2*l+1; m1++) {
|
||||
aa2 = m1 - l;
|
||||
for (int m2 = MAX(0,l-m1); m2 < MIN(2*l+1,3*l-m1+1); m2++) {
|
||||
bb2 = m2 - l;
|
||||
m = aa2 + bb2 + l;
|
||||
for (int m1 = 0; m1 < 2 * l + 1; m1++) {
|
||||
aa2 = m1 - l;
|
||||
for (int m2 = MAX(0, l - m1); m2 < MIN(2 * l + 1, 3 * l - m1 + 1); m2++) {
|
||||
bb2 = m2 - l;
|
||||
m = aa2 + bb2 + l;
|
||||
|
||||
sum = 0.0;
|
||||
for (int z = MAX(0, MAX(-aa2, bb2));
|
||||
z <= MIN(l, MIN(l - aa2, l + bb2)); z++) {
|
||||
ifac = z % 2 ? -1 : 1;
|
||||
sum += ifac /
|
||||
(factorial(z) *
|
||||
factorial(l - z) *
|
||||
factorial(l - aa2 - z) *
|
||||
factorial(l + bb2 - z) *
|
||||
factorial(aa2 + z) *
|
||||
factorial(-bb2 + z));
|
||||
}
|
||||
// clang-format off
|
||||
|
||||
cc2 = m - l;
|
||||
sfaccg = sqrt(factorial(l + aa2) *
|
||||
factorial(l - aa2) *
|
||||
factorial(l + bb2) *
|
||||
factorial(l - bb2) *
|
||||
factorial(l + cc2) *
|
||||
factorial(l - cc2) *
|
||||
(2*l + 1));
|
||||
|
||||
sfac1 = factorial(3*l + 1);
|
||||
sfac2 = factorial(l);
|
||||
dcg = sqrt(sfac2*sfac2*sfac2 / sfac1);
|
||||
|
||||
cglist[idxcg_count] = sum * dcg * sfaccg;
|
||||
idxcg_count++;
|
||||
sum = 0.0;
|
||||
for (int z = MAX(0, MAX(-aa2, bb2));
|
||||
z <= MIN(l, MIN(l - aa2, l + bb2)); z++) {
|
||||
ifac = z % 2 ? -1 : 1;
|
||||
sum += ifac /
|
||||
(factorial(z) *
|
||||
factorial(l - z) *
|
||||
factorial(l - aa2 - z) *
|
||||
factorial(l + bb2 - z) *
|
||||
factorial(aa2 + z) *
|
||||
factorial(-bb2 + z));
|
||||
}
|
||||
|
||||
cc2 = m - l;
|
||||
sfaccg = sqrt(factorial(l + aa2) *
|
||||
factorial(l - aa2) *
|
||||
factorial(l + bb2) *
|
||||
factorial(l - bb2) *
|
||||
factorial(l + cc2) *
|
||||
factorial(l - cc2) *
|
||||
(2*l + 1));
|
||||
// clang-format on
|
||||
|
||||
sfac1 = factorial(3 * l + 1);
|
||||
sfac2 = factorial(l);
|
||||
dcg = sqrt(sfac2 * sfac2 * sfac2 / sfac1);
|
||||
|
||||
cglist[idxcg_count] = sum * dcg * sfaccg;
|
||||
idxcg_count++;
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user