#include "fix_hyper_local.h"
#include <mpi.h>
#include <cmath>
#include <cstring>
#include "atom.h"
#include "update.h"
#include "force.h"
#include "pair.h"
#include "domain.h"
#include "comm.h"
#include "my_page.h"
#include "neighbor.h"
#include "neigh_request.h"
#include "neigh_list.h"
#include "math_extra.h"
#include "memory.h"
#include "error.h"
using namespace LAMMPS_NS;
using namespace FixConst;
#define DELTABOND 16384
#define DELTABIAS 16
#define COEFFINIT 1.0
#define FCCBONDS 12
#define BIG 1.0e20
enum{STRAIN,STRAINDOMAIN,BIASFLAG,BIASCOEFF};
enum{IGNORE,WARN,ERROR};
FixHyperLocal::FixHyperLocal(LAMMPS *lmp, int narg, char **arg) :
FixHyper(lmp, narg, arg), blist(NULL), biascoeff(NULL), numbond(NULL),
maxhalf(NULL), eligible(NULL), maxhalfstrain(NULL), old2now(NULL),
tagold(NULL), xold(NULL), maxstrain(NULL), maxstrain_domain(NULL),
biasflag(NULL), bias(NULL), cpage(NULL), clist(NULL), numcoeff(NULL)
{
if (atom->map_style == 0)
error->all(FLERR,"Fix hyper/local command requires atom map");
if (narg < 10) error->all(FLERR,"Illegal fix hyper/local command");
hyperflag = 2;
scalar_flag = 1;
vector_flag = 1;
size_vector = 21;
local_flag = 1;
size_local_rows = 0;
size_local_cols = 0;
local_freq = 1;
global_freq = 1;
extscalar = 0;
extvector = 0;
cutbond = force->numeric(FLERR,arg[3]);
qfactor = force->numeric(FLERR,arg[4]);
vmax = force->numeric(FLERR,arg[5]);
tequil = force->numeric(FLERR,arg[6]);
dcut = force->numeric(FLERR,arg[7]);
alpha_user = force->numeric(FLERR,arg[8]);
boost_target = force->numeric(FLERR,arg[9]);
if (cutbond < 0.0 || qfactor < 0.0 || vmax < 0.0 ||
tequil <= 0.0 || dcut <= 0.0 || alpha_user <= 0.0 || boost_target < 1.0)
error->all(FLERR,"Illegal fix hyper/local command");
invqfactorsq = 1.0 / (qfactor*qfactor);
cutbondsq = cutbond*cutbond;
dcutsq = dcut*dcut;
beta = 1.0 / (force->boltz * tequil);
checkghost = 0;
checkbias = 0;
int iarg = 10;
while (iarg < narg) {
if (strcmp(arg[iarg],"check/ghost") == 0) {
checkghost = 1;
iarg++;
} else if (strcmp(arg[iarg],"check/bias") == 0) {
if (iarg+3 > narg) error->all(FLERR,"Illegal fix hyper/local command");
checkbias = 1;
checkbias_every = force->inumeric(FLERR,arg[iarg+1]);
if (strcmp(arg[iarg+2],"error") == 0) checkbias_flag = ERROR;
else if (strcmp(arg[iarg+2],"warn") == 0) checkbias_flag = WARN;
else if (strcmp(arg[iarg+2],"ignore") == 0) checkbias_flag = IGNORE;
else error->all(FLERR,"Illegal fix hyper/local command");
iarg += 3;
} else error->all(FLERR,"Illegal fix hyper/local command");
}
maxbond = nblocal = 0;
blist = NULL;
biascoeff = NULL;
allbonds = 0;
maxatom = 0;
maxstrain = NULL;
maxstrain_domain = NULL;
biasflag = NULL;
maxlocal = nlocal_old = 0;
numbond = NULL;
maxhalf = NULL;
eligible = NULL;
maxhalfstrain = NULL;
maxall = nall_old = 0;
xold = NULL;
tagold = NULL;
old2now = NULL;
nbias = maxbias = 0;
bias = NULL;
maxcoeff = 0;
maxbondperatom = FCCBONDS;
numcoeff = NULL;
clist = NULL;
cpage = new MyPage<OneCoeff>;
cpage->init(maxbondperatom,1024*maxbondperatom,1);
comm_forward = 1;
comm_reverse = 2;
me = comm->me;
firstflag = 1;
sumbiascoeff = 0.0;
avebiascoeff = 0.0;
starttime = update->ntimestep;
nostrainyet = 1;
nnewbond = 0;
nevent = 0;
nevent_atom = 0;
mybias = 0.0;
}
FixHyperLocal::~FixHyperLocal()
{
memory->destroy(blist);
memory->destroy(biascoeff);
memory->destroy(maxstrain);
memory->destroy(maxstrain_domain);
memory->destroy(biasflag);
memory->destroy(numbond);
memory->destroy(maxhalf);
memory->destroy(eligible);
memory->destroy(maxhalfstrain);
memory->destroy(xold);
memory->destroy(tagold);
memory->destroy(old2now);
memory->destroy(bias);
memory->destroy(numcoeff);
memory->sfree(clist);
delete cpage;
}
int FixHyperLocal::setmask()
{
int mask = 0;
mask |= PRE_NEIGHBOR;
mask |= PRE_REVERSE;
mask |= MIN_PRE_NEIGHBOR;
mask |= THERMO_ENERGY;
return mask;
}
void FixHyperLocal::init_hyper()
{
ghost_toofar = 0;
checkbias_count = 0;
maxdriftsq = 0.0;
maxbondlen = 0.0;
avebiascoeff = 0.0;
minbiascoeff = BIG;
maxbiascoeff = 0.0;
nbias_running = 0;
nobias_running = 0;
negstrain_running = 0;
rmaxever = 0.0;
rmaxeverbig = 0.0;
nbondbuild = 0;
time_bondbuild = 0.0;
}
void FixHyperLocal::init()
{
if (force->newton_pair == 0)
error->all(FLERR,"Hyper local requires newton pair on");
if (atom->molecular && me == 0)
error->warning(FLERR,"Hyper local for molecular systems "
"requires care in defining hyperdynamic bonds");
if (firstflag) {
double cutghost;
if (force->pair)
cutghost = MAX(force->pair->cutforce+neighbor->skin,comm->cutghostuser);
else
cutghost = comm->cutghostuser;
if (cutghost < dcut)
error->all(FLERR,"Fix hyper/local domain cutoff exceeds ghost atom range - "
"use comm_modify cutoff command");
if (cutghost < dcut+cutbond/2.0 && me == 0)
error->warning(FLERR,"Fix hyper/local ghost atom range "
"may not allow for atom drift between events");
}
alpha = update->dt / alpha_user;
int irequest_full = neighbor->request(this,instance_me);
neighbor->requests[irequest_full]->id = 1;
neighbor->requests[irequest_full]->pair = 0;
neighbor->requests[irequest_full]->fix = 1;
neighbor->requests[irequest_full]->half = 0;
neighbor->requests[irequest_full]->full = 1;
neighbor->requests[irequest_full]->cut = 1;
neighbor->requests[irequest_full]->cutoff = dcut;
neighbor->requests[irequest_full]->occasional = 1;
int irequest_half = neighbor->request(this,instance_me);
neighbor->requests[irequest_half]->id = 2;
neighbor->requests[irequest_half]->pair = 0;
neighbor->requests[irequest_half]->fix = 1;
neighbor->requests[irequest_half]->occasional = 1;
}
void FixHyperLocal::init_list(int id, NeighList *ptr)
{
if (id == 1) listfull = ptr;
else if (id == 2) listhalf = ptr;
}
void FixHyperLocal::setup_pre_neighbor()
{
pre_neighbor();
}
void FixHyperLocal::setup_pre_reverse(int eflag, int vflag)
{
setupflag = 1;
pre_reverse(eflag,vflag);
setupflag = 0;
}
void FixHyperLocal::pre_neighbor()
{
int i,m,iold,jold,ilocal,jlocal;
for (i = 0; i < nall_old; i++) old2now[i] = -1;
for (m = 0; m < nblocal; m++) {
iold = blist[m].iold;
jold = blist[m].jold;
ilocal = old2now[iold];
jlocal = old2now[jold];
if (ilocal < 0) {
ilocal = atom->map(tagold[iold]);
ilocal = domain->closest_image(xold[iold],ilocal);
if (ilocal < 0)
error->one(FLERR,"Fix hyper/local bond atom not found");
old2now[iold] = ilocal;
}
if (jlocal < 0) {
jlocal = atom->map(tagold[jold]);
jlocal = domain->closest_image(xold[iold],jlocal); if (jlocal < 0)
error->one(FLERR,"Fix hyper/local bond atom not found");
old2now[jold] = jlocal;
}
blist[m].i = ilocal;
blist[m].j = jlocal;
}
for (iold = 0; iold < nall_old; iold++) {
if (old2now[iold] >= 0) continue;
if (tagold[iold] == 0) continue;
ilocal = atom->map(tagold[iold]);
old2now[iold] = ilocal;
if (ilocal < 0) ghost_toofar++;
}
}
void FixHyperLocal::pre_reverse(int , int )
{
int i,j,m,ii,jj,inum,jnum,iold,jold,ibond,nbond,ijhalf,ncount;
double delx,dely,delz;
double r,r0,estrain,emax,ebias,vbias,fbias,fbiasr;
double halfstrain,selfstrain;
int *ilist,*jlist,*numneigh,**firstneigh;
nostrainyet = 0;
int nlocal = atom->nlocal;
int nall = nlocal + atom->nghost;
if (maxatom < nall) {
memory->destroy(maxstrain);
memory->destroy(maxstrain_domain);
if (checkbias) memory->destroy(biasflag);
maxatom = atom->nmax;
memory->create(maxstrain,maxatom,"hyper/local:maxstrain");
memory->create(maxstrain_domain,maxatom,"hyper/local:maxstrain_domain");
if (checkbias) memory->create(biasflag,maxatom,"hyper/local:biasflag");
}
for (iold = 0; iold < nlocal_old; iold++) eligible[iold] = 1;
for (i = 0; i < nall; i++) maxstrain[i] = 0.0;
double **x = atom->x;
m = 0;
for (iold = 0; iold < nlocal_old; iold++) {
nbond = numbond[iold];
if (!nbond) {
eligible[iold] = 0;
continue;
}
halfstrain = 0.0;
for (ibond = 0; ibond < nbond; ibond++) {
i = blist[m].i;
j = blist[m].j;
delx = x[i][0] - x[j][0];
dely = x[i][1] - x[j][1];
delz = x[i][2] - x[j][2];
r = sqrt(delx*delx + dely*dely + delz*delz);
maxbondlen = MAX(r,maxbondlen);
r0 = blist[m].r0;
estrain = fabs(r-r0) / r0;
maxstrain[i] = MAX(maxstrain[i],estrain);
maxstrain[j] = MAX(maxstrain[j],estrain);
if (estrain > halfstrain) {
halfstrain = estrain;
ijhalf = m;
}
m++;
}
maxhalf[iold] = ijhalf;
maxhalfstrain[iold] = halfstrain;
}
commflag = STRAIN;
comm->reverse_comm_fix(this);
comm->forward_comm_fix(this);
for (i = 0; i < nall; i++) maxstrain_domain[i] = 0.0;
inum = listfull->inum;
ilist = listfull->ilist;
numneigh = listfull->numneigh;
firstneigh = listfull->firstneigh;
double rmax = rmaxever;
double rmaxbig = rmaxeverbig;
double *sublo = domain->sublo;
double *subhi = domain->subhi;
for (ii = 0; ii < inum; ii++) {
iold = ilist[ii];
if (eligible[iold] == 0) continue;
jlist = firstneigh[iold];
jnum = numneigh[iold];
i = old2now[iold];
emax = selfstrain = maxstrain[i];
ncount = 0;
for (jj = 0; jj < jnum; jj++) {
jold = jlist[jj];
j = old2now[jold];
if (j < 0) {
emax = MAX(emax,qfactor);
if (selfstrain == qfactor) ncount++;
continue;
}
emax = MAX(emax,maxstrain[j]);
if (selfstrain == maxstrain[j]) ncount++;
if (checkghost) {
if (j >= nlocal) {
if (x[j][0] < sublo[0]) rmaxbig = MAX(rmaxbig,sublo[0]-x[j][0]);
if (x[j][1] < sublo[1]) rmaxbig = MAX(rmaxbig,sublo[1]-x[j][1]);
if (x[j][2] < sublo[2]) rmaxbig = MAX(rmaxbig,sublo[2]-x[j][2]);
if (x[j][0] > subhi[0]) rmaxbig = MAX(rmaxbig,x[j][0]-subhi[0]);
if (x[j][1] > subhi[1]) rmaxbig = MAX(rmaxbig,x[j][1]-subhi[1]);
if (x[j][2] > subhi[2]) rmaxbig = MAX(rmaxbig,x[j][2]-subhi[2]);
if (maxstrain[j] < qfactor) {
if (x[j][0] < sublo[0]) rmax = MAX(rmax,sublo[0]-x[j][0]);
if (x[j][1] < sublo[1]) rmax = MAX(rmax,sublo[1]-x[j][1]);
if (x[j][2] < sublo[2]) rmax = MAX(rmax,sublo[2]-x[j][2]);
if (x[j][0] > subhi[0]) rmax = MAX(rmax,x[j][0]-subhi[0]);
if (x[j][1] > subhi[1]) rmax = MAX(rmax,x[j][1]-subhi[1]);
if (x[j][2] > subhi[2]) rmax = MAX(rmax,x[j][2]-subhi[2]);
}
}
}
}
if (maxhalfstrain[iold] < selfstrain) eligible[iold] = 0;
if (selfstrain < emax) eligible[iold] = 0;
else if (ncount > 1) {
eligible[iold] = 0;
emax = -emax;
}
maxstrain_domain[i] = emax;
}
commflag = STRAINDOMAIN;
comm->reverse_comm_fix(this);
comm->forward_comm_fix(this);
nbias = 0;
for (iold = 0; iold < nlocal_old; iold++) {
if (eligible[iold] == 0) continue;
j = blist[maxhalf[iold]].j;
if (maxstrain[j] != maxstrain_domain[j]) continue;
if (nbias == maxbias) {
maxbias += DELTABIAS;
memory->grow(bias,maxbias,"hyper/local:bias");
}
bias[nbias++] = maxhalf[iold];
}
double **f = atom->f;
int nobias = 0;
int negstrain = 0;
mybias = 0.0;
for (int ibias = 0; ibias < nbias; ibias++) {
m = bias[ibias];
i = blist[m].i;
j = blist[m].j;
if (maxstrain[i] >= qfactor) {
nobias++;
continue;
}
delx = x[i][0] - x[j][0];
dely = x[i][1] - x[j][1];
delz = x[i][2] - x[j][2];
r = sqrt(delx*delx + dely*dely + delz*delz);
r0 = blist[m].r0;
ebias = (r-r0) / r0;
vbias = biascoeff[m] * vmax * (1.0 - ebias*ebias*invqfactorsq);
fbias = biascoeff[m] * 2.0 * vmax * ebias * invqfactorsq;
fbiasr = fbias / r0 / r;
f[i][0] += delx*fbiasr;
f[i][1] += dely*fbiasr;
f[i][2] += delz*fbiasr;
f[j][0] -= delx*fbiasr;
f[j][1] -= dely*fbiasr;
f[j][2] -= delz*fbiasr;
if (ebias < 0.0) negstrain++;
mybias += vbias;
}
if (setupflag) return;
nbias_running += nbias;
nobias_running += nobias;
negstrain_running += negstrain;
double emaxi,emaxj,maxboost_domain,bc;
double mybiascoeff = 0.0;
for (m = 0; m < nblocal; m++) {
i = blist[m].i;
j = blist[m].j;
emaxi = maxstrain_domain[i];
emaxj = maxstrain_domain[j];
emax = MAX(emaxi,emaxj);
if (emax < qfactor) vbias = vmax * (1.0 - emax*emax*invqfactorsq);
else vbias = 0.0;
maxboost_domain = exp(beta * biascoeff[m]*vbias);
biascoeff[m] -= alpha * (maxboost_domain-boost_target) / boost_target;
bc = biascoeff[m];
mybiascoeff += bc;
minbiascoeff = MIN(minbiascoeff,bc);
maxbiascoeff = MAX(maxbiascoeff,bc);
}
MPI_Allreduce(&mybiascoeff,&sumbiascoeff,1,MPI_DOUBLE,MPI_SUM,world);
if (allbonds) avebiascoeff += sumbiascoeff/allbonds;
if (checkghost) {
double rmax2[2],rmax2all[2];
rmax2[0] = rmax;
rmax2[1] = rmaxbig;
MPI_Allreduce(&rmax2,&rmax2all,2,MPI_DOUBLE,MPI_MAX,world);
rmaxever = rmax2all[0];
rmaxeverbig = rmax2all[1];
}
if (checkbias && update->ntimestep % checkbias_every == 0) {
for (i = 0; i < nall; i++) biasflag[i] = 0;
tagint *tag = atom->tag;
for (int ibias = 0; ibias < nbias; ibias++) {
m = bias[ibias];
i = blist[m].i;
j = blist[m].j;
biasflag[i] = tag[j];
biasflag[j] = tag[i];
}
commflag = BIASFLAG;
comm->reverse_comm_fix(this);
comm->forward_comm_fix(this);
for (ii = 0; ii < inum; ii++) {
iold = ilist[ii];
i = old2now[iold];
if (biasflag[i] == 0) continue;
jlist = firstneigh[iold];
jnum = numneigh[iold];
for (jj = 0; jj < jnum; jj++) {
jold = jlist[jj];
j = old2now[jold];
if (j < 0) continue;
if (biasflag[j] && biasflag[j] != tag[i]) checkbias_count++;
}
}
if (checkbias_flag != IGNORE) {
int allcount;
MPI_Allreduce(&checkbias_count,&allcount,1,MPI_INT,MPI_SUM,world);
if (allcount) {
char str[128];
sprintf(str,"Fix hyper/local biased bonds too close: "
"cumulative atom count %d",allcount);
if (checkbias_flag == WARN) {
if (me == 0) error->warning(FLERR,str);
} else error->all(FLERR,str);
}
}
}
}
void FixHyperLocal::min_pre_neighbor()
{
pre_neighbor();
}
void FixHyperLocal::build_bond_list(int natom)
{
int i,j,ii,jj,m,n,iold,jold,ilocal,jlocal,inum,jnum,nbond;
tagint jtag;
double xtmp,ytmp,ztmp,delx,dely,delz,rsq,distsq,oldcoeff;
int *ilist,*jlist,*numneigh,**firstneigh;
double time1,time2;
time1 = MPI_Wtime();
if (natom) {
nevent++;
nevent_atom += natom;
}
double **x = atom->x;
for (i = 0; i < nall_old; i++) old2now[i] = -1;
for (m = 0; m < nblocal; m++) {
iold = blist[m].iold;
if (old2now[iold] < 0) {
ilocal = atom->map(tagold[iold]);
ilocal = domain->closest_image(xold[iold],ilocal);
if (ilocal < 0) error->one(FLERR,"Fix hyper/local bond atom not found");
old2now[iold] = ilocal;
distsq = MathExtra::distsq3(x[ilocal],xold[iold]);
maxdriftsq = MAX(distsq,maxdriftsq);
}
jold = blist[m].jold;
if (old2now[jold] < 0) {
jold = blist[m].jold;
jlocal = atom->map(tagold[jold]);
jlocal = domain->closest_image(xold[iold],jlocal); if (jlocal < 0) error->one(FLERR,"Fix hyper/local bond atom not found");
old2now[jold] = jlocal;
distsq = MathExtra::distsq3(x[jlocal],xold[jold]);
maxdriftsq = MAX(distsq,maxdriftsq);
}
}
tagint *tag = atom->tag;
int nlocal = atom->nlocal;
int nall = nlocal + atom->nghost;
if (maxcoeff < nall) {
memory->destroy(numcoeff);
memory->sfree(clist);
maxcoeff = atom->nmax;
memory->create(numcoeff,maxcoeff,"hyper/local:numcoeff");
clist = (OneCoeff **) memory->smalloc(maxcoeff*sizeof(OneCoeff *),
"hyper/local:clist");
}
while (1) {
if (firstflag) break;
for (i = 0; i < nall; i++) numcoeff[i] = 0;
for (i = 0; i < nall; i++) clist[i] = NULL;
cpage->reset();
for (m = 0; m < nblocal; m++) {
i = blist[m].i;
j = blist[m].j;
if (numcoeff[i] == 0) clist[i] = cpage->get(maxbondperatom);
if (numcoeff[j] == 0) clist[j] = cpage->get(maxbondperatom);
if (numcoeff[i] < maxbondperatom) {
clist[i][numcoeff[i]].biascoeff = biascoeff[m];
clist[i][numcoeff[i]].tag = tag[j];
}
numcoeff[i]++;
if (numcoeff[j] < maxbondperatom) {
clist[j][numcoeff[j]].biascoeff = biascoeff[m];
clist[j][numcoeff[j]].tag = tag[i];
}
numcoeff[j]++;
}
int mymax = 0;
for (i = 0; i < nall; i++) mymax = MAX(mymax,numcoeff[i]);
int maxcoeffall;
MPI_Allreduce(&mymax,&maxcoeffall,1,MPI_INT,MPI_MAX,world);
if (maxcoeffall > maxbondperatom) {
maxbondperatom = maxcoeffall;
cpage->init(maxbondperatom,1024*maxbondperatom,1);
continue;
}
commflag = BIASCOEFF;
comm->reverse_comm_fix_variable(this);
mymax = 0;
for (i = 0; i < nall; i++) mymax = MAX(mymax,numcoeff[i]);
MPI_Allreduce(&mymax,&maxcoeffall,1,MPI_INT,MPI_MAX,world);
if (maxcoeffall <= maxbondperatom) break;
maxbondperatom = maxcoeffall;
cpage->init(maxbondperatom,1024*maxbondperatom,1);
}
if (nlocal > maxlocal) {
memory->destroy(eligible);
memory->destroy(numbond);
memory->destroy(maxhalf);
memory->destroy(maxhalfstrain);
maxlocal = nlocal;
memory->create(eligible,maxlocal,"hyper/local:eligible");
memory->create(numbond,maxlocal,"hyper/local:numbond");
memory->create(maxhalf,maxlocal,"hyper/local:maxhalf");
memory->create(maxhalfstrain,maxlocal,"hyper/local:maxhalfstrain");
}
if (nall > maxall) {
memory->destroy(xold);
memory->destroy(tagold);
memory->destroy(old2now);
maxall = atom->nmax;
memory->create(xold,maxall,3,"hyper/local:xold");
memory->create(tagold,maxall,"hyper/local:tagold");
memory->create(old2now,maxall,"hyper/local:old2now");
}
nlocal_old = nlocal;
nall_old = nall;
memcpy(&xold[0][0],&x[0][0],3*nall*sizeof(double));
for (i = 0; i < nall; i++) tagold[i] = 0;
for (i = 0; i < nlocal; i++) numbond[i] = 0;
neighbor->build_one(listfull);
neighbor->build_one(listhalf);
if (listfull->inum != nlocal || listhalf->inum != nlocal)
error->one(FLERR,"Invalid neighbor list in fix hyper/local bond build");
inum = listfull->inum;
ilist = listfull->ilist;
numneigh = listfull->numneigh;
firstneigh = listfull->firstneigh;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
tagold[i] = tag[i];
jlist = firstneigh[i];
jnum = numneigh[i];
for (jj = 0; jj < jnum; jj++) {
j = jlist[jj];
j &= NEIGHMASK;
tagold[j] = tag[j];
}
}
int *mask = atom->mask;
inum = listhalf->inum;
ilist = listhalf->ilist;
numneigh = listhalf->numneigh;
firstneigh = listhalf->firstneigh;
nblocal = 0;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
xtmp = x[i][0];
ytmp = x[i][1];
ztmp = x[i][2];
tagold[i] = tag[i];
jlist = firstneigh[i];
jnum = numneigh[i];
nbond = 0;
for (jj = 0; jj < jnum; jj++) {
j = jlist[jj];
j &= NEIGHMASK;
jtag = tag[j];
tagold[j] = jtag;
if (!(mask[i] & groupbit) && !(mask[j] & groupbit)) continue;
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
delz = ztmp - x[j][2];
rsq = delx*delx + dely*dely + delz*delz;
if (rsq < cutbondsq) {
nbond++;
if (nblocal == maxbond) grow_bond();
blist[nblocal].i = i;
blist[nblocal].j = j;
blist[nblocal].iold = i;
blist[nblocal].jold = j;
blist[nblocal].r0 = sqrt(rsq);
if (firstflag) oldcoeff = 0.0;
else {
oldcoeff = 0.0;
jtag = tag[j];
n = numcoeff[i];
for (m = 0; m < n; m++) {
if (clist[i][m].tag == jtag) {
oldcoeff = clist[i][m].biascoeff;
break;
}
}
}
if (oldcoeff > 0.0) biascoeff[nblocal] = oldcoeff;
else {
biascoeff[nblocal] = COEFFINIT;
nnewbond++;
}
nblocal++;
}
}
numbond[i] = nbond;
}
size_local_rows = nblocal;
bigint bondcount = nblocal;
MPI_Allreduce(&bondcount,&allbonds,1,MPI_LMP_BIGINT,MPI_SUM,world);
time2 = MPI_Wtime();
if (firstflag) nnewbond = 0;
else {
time_bondbuild += time2-time1;
nbondbuild++;
}
firstflag = 0;
}
int FixHyperLocal::pack_forward_comm(int n, int *list, double *buf,
int , int * )
{
int i,j,m;
m = 0;
if (commflag == STRAIN) {
for (i = 0; i < n; i++) {
j = list[i];
buf[m++] = maxstrain[j];
}
} else if (commflag == STRAINDOMAIN) {
for (i = 0; i < n; i++) {
j = list[i];
buf[m++] = maxstrain_domain[j];
}
} else if (commflag == BIASFLAG) {
for (i = 0; i < n; i++) {
j = list[i];
buf[m++] = ubuf(biasflag[j]).d;
}
}
return m;
}
void FixHyperLocal::unpack_forward_comm(int n, int first, double *buf)
{
int i,m,last;
m = 0;
last = first + n;
if (commflag == STRAIN) {
for (i = first; i < last; i++) {
maxstrain[i] = buf[m++];
}
} else if (commflag == STRAINDOMAIN) {
for (i = first; i < last; i++) {
maxstrain_domain[i] = buf[m++];
}
} else if (commflag == BIASFLAG) {
for (i = first; i < last; i++) {
biasflag[i] = (tagint) ubuf(buf[m++]).i;
}
}
}
int FixHyperLocal::pack_reverse_comm(int n, int first, double *buf)
{
int i,j,m,last;
m = 0;
last = first + n;
if (commflag == STRAIN) {
int nonzero = 0;
m++; for (i = first; i < last; i++) {
if (maxstrain[i] == 0.0) continue;
nonzero++;
buf[m++] = ubuf(i-first).d; buf[m++] = maxstrain[i]; }
buf[0] = ubuf(nonzero).d;
} else if (commflag == STRAINDOMAIN) {
int nonzero = 0;
m++; for (i = first; i < last; i++) {
if (maxstrain_domain[i] == 0.0) continue;
nonzero++;
buf[m++] = ubuf(i-first).d; buf[m++] = maxstrain_domain[i]; }
buf[0] = ubuf(nonzero).d;
} else if (commflag == BIASFLAG) {
for (i = first; i < last; i++) {
buf[m++] = ubuf(biasflag[i]).d;
}
} else if (commflag == BIASCOEFF) {
int ncoeff;
int nonzero = 0;
m++; for (i = first; i < last; i++) {
if (numcoeff[i] == 0) continue;
nonzero++;
ncoeff = numcoeff[i];
buf[m++] = ubuf(i-first).d; buf[m++] = ubuf(ncoeff).d; for (j = 0; j < ncoeff; j++) {
buf[m++] = clist[i][j].biascoeff;
buf[m++] = ubuf(clist[i][j].tag).d;
}
}
buf[0] = ubuf(nonzero).d;
}
return m;
}
int FixHyperLocal::pack_reverse_comm_size(int n, int first)
{
int last = first + n;
int m = 1;
for (int i = first; i < last; i++) {
if (numcoeff[i]) m += 2 + 2*numcoeff[i];
}
return m;
}
void FixHyperLocal::unpack_reverse_comm(int n, int *list, double *buf)
{
int i,j,k,m;
if (n == 0) return;
m = 0;
if (commflag == STRAIN) {
int offset;
int nonzero = (int) ubuf(buf[m++]).i;
for (int iatom = 0; iatom < nonzero; iatom++) {
offset = (int) ubuf(buf[m++]).i; j = list[offset];
maxstrain[j] = MAX(maxstrain[j],buf[m]);
m++;
}
} else if (commflag == STRAINDOMAIN) {
int offset;
int nonzero = (int) ubuf(buf[m++]).i; for (int iatom = 0; iatom < nonzero; iatom++) {
offset = (int) ubuf(buf[m++]).i; j = list[offset];
maxstrain_domain[j] = MAX(maxstrain_domain[j],buf[m]);
m++;
}
} else if (commflag == BIASFLAG) {
for (i = 0; i < n; i++) {
j = list[i];
biasflag[j] = (tagint) ubuf(buf[m++]).i;
}
} else if (commflag == BIASCOEFF) {
int offset,ncoeff;
int nonzero = (int) ubuf(buf[m++]).i; for (int iatom = 0; iatom < nonzero; iatom++) {
offset = (int) ubuf(buf[m++]).i; j = list[offset];
ncoeff = (int) ubuf(buf[m++]).i; for (k = 0; k < ncoeff; k++) {
if (numcoeff[j] == 0) clist[j] = cpage->get(maxbondperatom);
if (numcoeff[j] < maxbondperatom) {
clist[j][numcoeff[j]].biascoeff = buf[m++];
clist[j][numcoeff[j]].tag = (tagint) ubuf(buf[m++]).i;
} else m += 2;
numcoeff[j]++;
}
}
}
}
void FixHyperLocal::grow_bond()
{
if (maxbond + DELTABOND > MAXSMALLINT)
error->one(FLERR,"Fix hyper/local bond count is too big");
maxbond += DELTABOND;
blist = (OneBond *)
memory->srealloc(blist,maxbond*sizeof(OneBond),"hyper/local:blist");
memory->grow(biascoeff,maxbond,"hyper/local:biascoeff");
vector_local = biascoeff;
}
double FixHyperLocal::compute_scalar()
{
double allbias;
MPI_Allreduce(&mybias,&allbias,1,MPI_DOUBLE,MPI_SUM,world);
return allbias;
}
double FixHyperLocal::compute_vector(int i)
{
if (i == 0) {
int nbiasall;
MPI_Allreduce(&nbias,&nbiasall,1,MPI_INT,MPI_SUM,world);
return (double) nbiasall;
}
if (i == 1) {
if (nostrainyet) return 0.0;
int nlocal = atom->nlocal;
double emax = 0.0;
for (int j = 0; j < nlocal; j++)
emax = MAX(emax,maxstrain[j]);
double eall;
MPI_Allreduce(&emax,&eall,1,MPI_DOUBLE,MPI_MAX,world);
return eall;
}
if (i == 2) {
if (allbonds) return sumbiascoeff/allbonds;
return 1.0;
}
if (i == 3) {
return 2.0*allbonds/atom->natoms;
}
if (i == 4) {
bigint allneigh,thisneigh;
thisneigh = listfull->ipage->ndatum;
MPI_Allreduce(&thisneigh,&allneigh,1,MPI_LMP_BIGINT,MPI_SUM,world);
const double natoms = atom->natoms;
const double neighsperatom = static_cast<double>(allneigh)/natoms;
const double bondsperatom = static_cast<double>(allbonds)/natoms;
return neighsperatom * bondsperatom;
}
if (i == 5) {
double allbondlen;
MPI_Allreduce(&maxbondlen,&allbondlen,1,MPI_DOUBLE,MPI_MAX,world);
return allbondlen;
}
if (i == 6) {
if (update->ntimestep == update->firststep) return 0.0;
int allbias_running;
MPI_Allreduce(&nbias_running,&allbias_running,1,MPI_INT,MPI_SUM,world);
return 1.0*allbias_running / (update->ntimestep - update->firststep);
}
if (i == 7) {
int allbias_running,allnobias_running;
MPI_Allreduce(&nbias_running,&allbias_running,1,MPI_INT,MPI_SUM,world);
MPI_Allreduce(&nobias_running,&allnobias_running,1,MPI_INT,MPI_SUM,world);
if (allbias_running) return 1.0*allnobias_running / allbias_running;
return 0.0;
}
if (i == 8) {
int allbias_running,allnegstrain_running;
MPI_Allreduce(&nbias_running,&allbias_running,1,MPI_INT,MPI_SUM,world);
MPI_Allreduce(&negstrain_running,&allnegstrain_running,1,MPI_INT,
MPI_SUM,world);
if (allbias_running) return 1.0*allnegstrain_running / allbias_running;
return 0.0;
}
if (i == 9) {
if (update->ntimestep == update->firststep) return 0.0;
return avebiascoeff / (update->ntimestep - update->firststep);
}
if (i == 10) {
double biascoeff;
MPI_Allreduce(&minbiascoeff,&biascoeff,1,MPI_DOUBLE,MPI_MIN,world);
return biascoeff;
}
if (i == 11) {
double biascoeff;
MPI_Allreduce(&maxbiascoeff,&biascoeff,1,MPI_DOUBLE,MPI_MAX,world);
return biascoeff;
}
if (i == 12) {
double alldriftsq;
MPI_Allreduce(&maxdriftsq,&alldriftsq,1,MPI_DOUBLE,MPI_MAX,world);
return (double) sqrt(alldriftsq);
}
if (i == 13) return rmaxever;
if (i == 14) return rmaxeverbig;
if (i == 15) {
int allghost_toofar;
MPI_Allreduce(&ghost_toofar,&allghost_toofar,1,MPI_INT,MPI_SUM,world);
return 1.0*allghost_toofar;
}
if (i == 16) {
int allclose;
MPI_Allreduce(&checkbias_count,&allclose,1,MPI_INT,MPI_SUM,world);
return 1.0*allclose;
}
if (i == 17) {
return boost_target * update->dt * (update->ntimestep - starttime);
}
if (i == 18) return (double) nevent;
if (i == 19) return (double) nevent_atom;
if (i == 20) {
int allnewbond;
MPI_Allreduce(&nnewbond,&allnewbond,1,MPI_INT,MPI_SUM,world);
return (double) allnewbond;
}
return 0.0;
}
double FixHyperLocal::query(int i)
{
if (i == 1) return compute_vector(17); if (i == 2) return compute_vector(18); if (i == 3) return compute_vector(19); if (i == 4) return compute_vector(3); if (i == 5) return compute_vector(12); if (i == 6) return compute_vector(5); if (i == 7) return compute_vector(7); if (i == 8) return compute_vector(8);
if (i == 9) return compute_vector(20); if (i == 10) return 1.0*maxbondperatom; if (i == 11) return compute_vector(6); if (i == 12) return compute_vector(9); if (i == 13) return compute_vector(10); if (i == 14) return compute_vector(11); if (i == 15) return compute_vector(4); if (i == 16) return compute_vector(2); if (i == 17) return time_bondbuild; if (i == 18) return rmaxever; if (i == 19) return rmaxeverbig; if (i == 20) return compute_vector(15); if (i == 21) return compute_vector(16);
error->all(FLERR,"Invalid query to fix hyper/local");
return 0.0;
}
double FixHyperLocal::memory_usage()
{
double bytes = maxbond * sizeof(OneBond); bytes = maxbond * sizeof(double); bytes += 3*maxlocal * sizeof(int); bytes += maxlocal * sizeof(double); bytes += maxall * sizeof(int); bytes += maxall * sizeof(tagint); bytes += 3*maxall * sizeof(double); bytes += 2*maxall * sizeof(double); if (checkbias) bytes += maxall * sizeof(tagint); bytes += maxcoeff * sizeof(int); bytes += maxcoeff * sizeof(OneCoeff *); bytes += maxlocal*maxbondperatom * sizeof(OneCoeff); return bytes;
}