#include "special.h"
#include <mpi.h>
#include "atom.h"
#include "atom_vec.h"
#include "force.h"
#include "comm.h"
#include "modify.h"
#include "fix.h"
#include "accelerator_kokkos.h"
#include "atom_masks.h"
#include "memory.h"
using namespace LAMMPS_NS;
#define RVOUS 1
Special::Special(LAMMPS *lmp) : Pointers(lmp)
{
MPI_Comm_rank(world,&me);
MPI_Comm_size(world,&nprocs);
onetwo = onethree = onefour = NULL;
}
Special::~Special()
{
memory->destroy(onetwo);
memory->destroy(onethree);
memory->destroy(onefour);
}
void Special::build()
{
MPI_Barrier(world);
double time1 = MPI_Wtime();
if (me == 0 && screen) {
const double * const special_lj = force->special_lj;
const double * const special_coul = force->special_coul;
fprintf(screen,"Finding 1-2 1-3 1-4 neighbors ...\n"
" special bond factors lj: %-10g %-10g %-10g\n"
" special bond factors coul: %-10g %-10g %-10g\n",
special_lj[1],special_lj[2],special_lj[3],
special_coul[1],special_coul[2],special_coul[3]);
}
int **nspecial = atom->nspecial;
int nlocal = atom->nlocal;
for (int i = 0; i < nlocal; i++) {
nspecial[i][0] = 0;
nspecial[i][1] = 0;
nspecial[i][2] = 0;
}
atom_owners();
if (force->newton_bond) onetwo_build_newton();
else onetwo_build_newton_off();
if (me == 0) {
if (screen) fprintf(screen," %d = max # of 1-2 neighbors\n",maxall);
if (logfile) fprintf(logfile," %d = max # of 1-2 neighbors\n",maxall);
}
if (force->special_lj[2] == 1.0 && force->special_coul[2] == 1.0 &&
force->special_lj[3] == 1.0 && force->special_coul[3] == 1.0) {
dedup();
combine();
fix_alteration();
memory->destroy(procowner);
memory->destroy(atomIDs);
timer_output(time1);
return;
}
onethree_build();
if (me == 0) {
if (screen) fprintf(screen," %d = max # of 1-3 neighbors\n",maxall);
if (logfile) fprintf(logfile," %d = max # of 1-3 neighbors\n",maxall);
}
if (force->special_lj[3] == 1.0 && force->special_coul[3] == 1.0) {
dedup();
if (force->special_angle) angle_trim();
combine();
fix_alteration();
memory->destroy(procowner);
memory->destroy(atomIDs);
timer_output(time1);
return;
}
onefour_build();
if (me == 0) {
if (screen) fprintf(screen," %d = max # of 1-4 neighbors\n",maxall);
if (logfile) fprintf(logfile," %d = max # of 1-4 neighbors\n",maxall);
}
dedup();
if (force->special_angle) angle_trim();
if (force->special_dihedral) dihedral_trim();
combine();
fix_alteration();
memory->destroy(procowner);
memory->destroy(atomIDs);
timer_output(time1);
}
void Special::atom_owners()
{
tagint *tag = atom->tag;
int nlocal = atom->nlocal;
int *proclist;
memory->create(proclist,nlocal,"special:proclist");
IDRvous *idbuf = (IDRvous *)
memory->smalloc((bigint) nlocal*sizeof(IDRvous),"special:idbuf");
for (int i = 0; i < nlocal; i++) {
proclist[i] = tag[i] % nprocs;
idbuf[i].me = me;
idbuf[i].atomID = tag[i];
}
char *buf;
comm->rendezvous(RVOUS,nlocal,(char *) idbuf,sizeof(IDRvous),0,proclist,
rendezvous_ids,0,buf,0,(void *) this);
memory->destroy(proclist);
memory->sfree(idbuf);
}
void Special::onetwo_build_newton()
{
int i,j,m;
tagint *tag = atom->tag;
int *num_bond = atom->num_bond;
tagint **bond_atom = atom->bond_atom;
int **nspecial = atom->nspecial;
int nlocal = atom->nlocal;
int nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_bond[i]; j++) {
m = atom->map(bond_atom[i][j]);
if (m < 0 || m >= nlocal) nsend++;
}
}
int *proclist;
memory->create(proclist,nsend,"special:proclist");
PairRvous *inbuf = (PairRvous *)
memory->smalloc((bigint) nsend*sizeof(PairRvous),"special:inbuf");
nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_bond[i]; j++) {
m = atom->map(bond_atom[i][j]);
if (m >= 0 && m < nlocal) continue;
proclist[nsend] = bond_atom[i][j] % nprocs;
inbuf[nsend].atomID = bond_atom[i][j];
inbuf[nsend].partnerID = tag[i];
nsend++;
}
}
char *buf;
int nreturn = comm->rendezvous(RVOUS,nsend,(char *) inbuf,sizeof(PairRvous),
0,proclist,
rendezvous_pairs,0,buf,sizeof(PairRvous),
(void *) this);
PairRvous *outbuf = (PairRvous *) buf;
memory->destroy(proclist);
memory->sfree(inbuf);
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_bond[i]; j++) {
nspecial[i][0]++;
m = atom->map(bond_atom[i][j]);
if (m >= 0 && m < nlocal) nspecial[m][0]++;
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
nspecial[i][0]++;
}
int max = 0;
for (i = 0; i < nlocal; i++)
max = MAX(max,nspecial[i][0]);
MPI_Allreduce(&max,&maxall,1,MPI_INT,MPI_MAX,world);
memory->create(onetwo,nlocal,maxall,"special:onetwo");
for (i = 0; i < nlocal; i++) nspecial[i][0] = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_bond[i]; j++) {
onetwo[i][nspecial[i][0]++] = bond_atom[i][j];
m = atom->map(bond_atom[i][j]);
if (m >= 0 && m < nlocal) onetwo[m][nspecial[m][0]++] = tag[i];
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
onetwo[i][nspecial[i][0]++] = outbuf[m].partnerID;
}
memory->sfree(outbuf);
}
void Special::onetwo_build_newton_off()
{
int i,j;
int *num_bond = atom->num_bond;
tagint **bond_atom = atom->bond_atom;
int **nspecial = atom->nspecial;
int nlocal = atom->nlocal;
int max = 0;
for (i = 0; i < nlocal; i++)
max = MAX(max,num_bond[i]);
MPI_Allreduce(&max,&maxall,1,MPI_INT,MPI_MAX,world);
memory->create(onetwo,nlocal,maxall,"special:onetwo");
for (i = 0; i < nlocal; i++) {
nspecial[i][0] = num_bond[i];
for (j = 0; j < num_bond[i]; j++)
onetwo[i][j] = bond_atom[i][j];
}
}
void Special::onethree_build()
{
int i,j,k,m,proc;
int **nspecial = atom->nspecial;
int nlocal = atom->nlocal;
int nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][0]; j++) {
m = atom->map(onetwo[i][j]);
if (m < 0 || m >= nlocal) nsend += nspecial[i][0]-1;
}
}
int *proclist;
memory->create(proclist,nsend,"special:proclist");
PairRvous *inbuf = (PairRvous *)
memory->smalloc((bigint) nsend*sizeof(PairRvous),"special:inbuf");
nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][0]; j++) {
m = atom->map(onetwo[i][j]);
if (m >= 0 && m < nlocal) continue;
proc = onetwo[i][j] % nprocs;
for (k = 0; k < nspecial[i][0]; k++) {
if (j == k) continue;
proclist[nsend] = proc;
inbuf[nsend].atomID = onetwo[i][j];
inbuf[nsend].partnerID = onetwo[i][k];
nsend++;
}
}
}
char *buf;
int nreturn = comm->rendezvous(RVOUS,nsend,(char *) inbuf,sizeof(PairRvous),
0,proclist,
rendezvous_pairs,0,buf,sizeof(PairRvous),
(void *) this);
PairRvous *outbuf = (PairRvous *) buf;
memory->destroy(proclist);
memory->sfree(inbuf);
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][0]; j++) {
m = atom->map(onetwo[i][j]);
if (m >= 0 && m < nlocal) nspecial[m][1] += nspecial[i][0]-1;
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
nspecial[i][1]++;
}
int max = 0;
for (i = 0; i < nlocal; i++)
max = MAX(max,nspecial[i][1]);
MPI_Allreduce(&max,&maxall,1,MPI_INT,MPI_MAX,world);
memory->create(onethree,nlocal,maxall,"special:onethree");
for (i = 0; i < nlocal; i++) nspecial[i][1] = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][0]; j++) {
m = atom->map(onetwo[i][j]);
if (m < 0 || m >= nlocal) continue;
for (k = 0; k < nspecial[i][0]; k++) {
if (j == k) continue;
onethree[m][nspecial[m][1]++] = onetwo[i][k];
}
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
onethree[i][nspecial[i][1]++] = outbuf[m].partnerID;
}
memory->sfree(outbuf);
}
void Special::onefour_build()
{
int i,j,k,m,proc;
int **nspecial = atom->nspecial;
int nlocal = atom->nlocal;
int nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][1]; j++) {
m = atom->map(onethree[i][j]);
if (m < 0 || m >= nlocal) nsend += nspecial[i][0];
}
}
int *proclist;
memory->create(proclist,nsend,"special:proclist");
PairRvous *inbuf = (PairRvous *)
memory->smalloc((bigint) nsend*sizeof(PairRvous),"special:inbuf");
nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][1]; j++) {
m = atom->map(onethree[i][j]);
if (m >= 0 && m < nlocal) continue;
proc = onethree[i][j] % nprocs;
for (k = 0; k < nspecial[i][0]; k++) {
proclist[nsend] = proc;
inbuf[nsend].atomID = onethree[i][j];
inbuf[nsend].partnerID = onetwo[i][k];
nsend++;
}
}
}
char *buf;
int nreturn = comm->rendezvous(RVOUS,nsend,(char *) inbuf,sizeof(PairRvous),
0,proclist,
rendezvous_pairs,0,buf,sizeof(PairRvous),
(void *) this);
PairRvous *outbuf = (PairRvous *) buf;
memory->destroy(proclist);
memory->sfree(inbuf);
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][1]; j++) {
m = atom->map(onethree[i][j]);
if (m >= 0 && m < nlocal) nspecial[m][2] += nspecial[i][0];
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
nspecial[i][2]++;
}
int max = 0;
for (i = 0; i < nlocal; i++)
max = MAX(max,nspecial[i][2]);
MPI_Allreduce(&max,&maxall,1,MPI_INT,MPI_MAX,world);
memory->create(onefour,nlocal,maxall,"special:onefour");
for (i = 0; i < nlocal; i++) nspecial[i][2] = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < nspecial[i][1]; j++) {
m = atom->map(onethree[i][j]);
if (m < 0 || m >= nlocal) continue;
for (k = 0; k < nspecial[i][0]; k++) {
onefour[m][nspecial[m][2]++] = onetwo[i][k];
}
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
onefour[i][nspecial[i][2]++] = outbuf[m].partnerID;
}
memory->sfree(outbuf);
}
void Special::dedup()
{
int i,j;
tagint m;
atom->map_clear();
int **nspecial = atom->nspecial;
tagint *tag = atom->tag;
int nlocal = atom->nlocal;
int unique;
for (i = 0; i < nlocal; i++) {
unique = 0;
atom->map_one(tag[i],0);
for (j = 0; j < nspecial[i][0]; j++) {
m = onetwo[i][j];
if (atom->map(m) < 0) {
onetwo[i][unique++] = m;
atom->map_one(m,0);
}
}
nspecial[i][0] = unique;
atom->map_one(tag[i],-1);
for (j = 0; j < unique; j++) atom->map_one(onetwo[i][j],-1);
}
for (i = 0; i < nlocal; i++) {
unique = 0;
atom->map_one(tag[i],0);
for (j = 0; j < nspecial[i][1]; j++) {
m = onethree[i][j];
if (atom->map(m) < 0) {
onethree[i][unique++] = m;
atom->map_one(m,0);
}
}
nspecial[i][1] = unique;
atom->map_one(tag[i],-1);
for (j = 0; j < unique; j++) atom->map_one(onethree[i][j],-1);
}
for (i = 0; i < nlocal; i++) {
unique = 0;
atom->map_one(tag[i],0);
for (j = 0; j < nspecial[i][2]; j++) {
m = onefour[i][j];
if (atom->map(m) < 0) {
onefour[i][unique++] = m;
atom->map_one(m,0);
}
}
nspecial[i][2] = unique;
atom->map_one(tag[i],-1);
for (j = 0; j < unique; j++) atom->map_one(onefour[i][j],-1);
}
atom->map_init(0);
atom->nghost = 0;
atom->map_set();
}
void Special::combine()
{
int i,j;
tagint m;
int me;
MPI_Comm_rank(world,&me);
int **nspecial = atom->nspecial;
tagint *tag = atom->tag;
int nlocal = atom->nlocal;
atom->map_clear();
int unique;
int maxspecial = 0;
for (i = 0; i < nlocal; i++) {
unique = 0;
atom->map_one(tag[i],0);
for (j = 0; j < nspecial[i][0]; j++) {
m = onetwo[i][j];
if (atom->map(m) < 0) {
unique++;
atom->map_one(m,0);
}
}
for (j = 0; j < nspecial[i][1]; j++) {
m = onethree[i][j];
if (atom->map(m) < 0) {
unique++;
atom->map_one(m,0);
}
}
for (j = 0; j < nspecial[i][2]; j++) {
m = onefour[i][j];
if (atom->map(m) < 0) {
unique++;
atom->map_one(m,0);
}
}
maxspecial = MAX(maxspecial,unique);
atom->map_one(tag[i],-1);
for (j = 0; j < nspecial[i][0]; j++) atom->map_one(onetwo[i][j],-1);
for (j = 0; j < nspecial[i][1]; j++) atom->map_one(onethree[i][j],-1);
for (j = 0; j < nspecial[i][2]; j++) atom->map_one(onefour[i][j],-1);
}
maxspecial = MAX(atom->maxspecial,maxspecial);
MPI_Allreduce(&maxspecial,&atom->maxspecial,1,MPI_INT,MPI_MAX,world);
atom->maxspecial += force->special_extra;
force->special_extra = 0;
if (me == 0) {
if (screen)
fprintf(screen," %d = max # of special neighbors\n",atom->maxspecial);
if (logfile)
fprintf(logfile," %d = max # of special neighbors\n",atom->maxspecial);
}
if (lmp->kokkos) {
AtomKokkos* atomKK = (AtomKokkos*) atom;
atomKK->modified(Host,SPECIAL_MASK);
atomKK->sync(Device,SPECIAL_MASK);
MemoryKokkos* memoryKK = (MemoryKokkos*) memory;
memoryKK->grow_kokkos(atomKK->k_special,atom->special,
atom->nmax,atom->maxspecial,"atom:special");
atomKK->modified(Device,SPECIAL_MASK);
atomKK->sync(Host,SPECIAL_MASK);
} else {
memory->destroy(atom->special);
memory->create(atom->special,atom->nmax,atom->maxspecial,"atom:special");
}
atom->avec->grow_reset();
tagint **special = atom->special;
for (i = 0; i < nlocal; i++) {
unique = 0;
atom->map_one(tag[i],0);
for (j = 0; j < nspecial[i][0]; j++) {
m = onetwo[i][j];
if (atom->map(m) < 0) {
special[i][unique++] = m;
atom->map_one(m,0);
}
}
nspecial[i][0] = unique;
for (j = 0; j < nspecial[i][1]; j++) {
m = onethree[i][j];
if (atom->map(m) < 0) {
special[i][unique++] = m;
atom->map_one(m,0);
}
}
nspecial[i][1] = unique;
for (j = 0; j < nspecial[i][2]; j++) {
m = onefour[i][j];
if (atom->map(m) < 0) {
special[i][unique++] = m;
atom->map_one(m,0);
}
}
nspecial[i][2] = unique;
atom->map_one(tag[i],-1);
for (j = 0; j < nspecial[i][2]; j++) atom->map_one(special[i][j],-1);
}
atom->map_init(0);
atom->nghost = 0;
atom->map_set();
}
void Special::angle_trim()
{
int i,j,k,m;;
int *num_angle = atom->num_angle;
int *num_dihedral = atom->num_dihedral;
tagint **angle_atom1 = atom->angle_atom1;
tagint **angle_atom2 = atom->angle_atom2;
tagint **angle_atom3 = atom->angle_atom3;
tagint **dihedral_atom1 = atom->dihedral_atom1;
tagint **dihedral_atom2 = atom->dihedral_atom2;
tagint **dihedral_atom3 = atom->dihedral_atom3;
tagint **dihedral_atom4 = atom->dihedral_atom4;
int **nspecial = atom->nspecial;
tagint *tag = atom->tag;
int nlocal = atom->nlocal;
double onethreecount = 0.0;
for (i = 0; i < nlocal; i++) onethreecount += nspecial[i][1];
double allcount;
MPI_Allreduce(&onethreecount,&allcount,1,MPI_DOUBLE,MPI_SUM,world);
if (me == 0) {
if (screen)
fprintf(screen,
" %g = # of 1-3 neighbors before angle trim\n",allcount);
if (logfile)
fprintf(logfile,
" %g = # of 1-3 neighbors before angle trim\n",allcount);
}
if ((num_angle && atom->nangles) || (num_dihedral && atom->ndihedrals)) {
int nsend = 0;
for (i = 0; i < nlocal; i++) {
if (num_angle) {
for (j = 0; j < num_angle[i]; j++) {
if (tag[i] != angle_atom2[i][j]) continue;
m = atom->map(angle_atom1[i][j]);
if (m < 0 || m >= nlocal) nsend++;
m = atom->map(angle_atom3[i][j]);
if (m < 0 || m >= nlocal) nsend++;
}
}
if (num_dihedral) {
for (j = 0; j < num_dihedral[i]; j++) {
if (tag[i] != dihedral_atom2[i][j]) continue;
m = atom->map(dihedral_atom1[i][j]);
if (m < 0 || m >= nlocal) nsend++;
m = atom->map(dihedral_atom3[i][j]);
if (m < 0 || m >= nlocal) nsend++;
m = atom->map(dihedral_atom4[i][j]);
if (m < 0 || m >= nlocal) nsend++;
}
}
}
int *proclist;
memory->create(proclist,nsend,"special:proclist");
PairRvous *inbuf = (PairRvous *)
memory->smalloc((bigint) nsend*sizeof(PairRvous),"special:inbuf");
nsend = 0;
for (i = 0; i < nlocal; i++) {
if (num_angle) {
for (j = 0; j < num_angle[i]; j++) {
if (tag[i] != angle_atom2[i][j]) continue;
m = atom->map(angle_atom1[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = angle_atom1[i][j] % nprocs;
inbuf[nsend].atomID = angle_atom1[i][j];
inbuf[nsend].partnerID = angle_atom3[i][j];
nsend++;
}
m = atom->map(angle_atom3[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = angle_atom3[i][j] % nprocs;
inbuf[nsend].atomID = angle_atom3[i][j];
inbuf[nsend].partnerID = angle_atom1[i][j];
nsend++;
}
}
}
if (num_dihedral) {
for (j = 0; j < num_dihedral[i]; j++) {
if (tag[i] != dihedral_atom2[i][j]) continue;
m = atom->map(dihedral_atom1[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = dihedral_atom1[i][j] % nprocs;
inbuf[nsend].atomID = dihedral_atom1[i][j];
inbuf[nsend].partnerID = dihedral_atom3[i][j];
nsend++;
}
m = atom->map(dihedral_atom3[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = dihedral_atom3[i][j] % nprocs;
inbuf[nsend].atomID = dihedral_atom3[i][j];
inbuf[nsend].partnerID = dihedral_atom1[i][j];
nsend++;
}
m = atom->map(dihedral_atom4[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = dihedral_atom4[i][j] % nprocs;
inbuf[nsend].atomID = dihedral_atom4[i][j];
inbuf[nsend].partnerID = dihedral_atom2[i][j];
nsend++;
}
}
}
}
char *buf;
int nreturn = comm->rendezvous(RVOUS,nsend,(char *) inbuf,sizeof(PairRvous),
0,proclist,
rendezvous_pairs,0,buf,sizeof(PairRvous),
(void *) this);
PairRvous *outbuf = (PairRvous *) buf;
memory->destroy(proclist);
memory->sfree(inbuf);
int max = 0;
for (i = 0; i < nlocal; i++)
max = MAX(max,nspecial[i][1]);
MPI_Allreduce(&max,&maxall,1,MPI_INT,MPI_MAX,world);
int **flag;
memory->create(flag,nlocal,maxall,"special:flag");
for (i = 0; i < nlocal; i++)
for (j = 0; j < nspecial[i][1]; j++)
flag[i][j] = 0;
for (i = 0; i < nlocal; i++) {
if (num_angle) {
for (j = 0; j < num_angle[i]; j++) {
if (tag[i] != angle_atom2[i][j]) continue;
m = atom->map(angle_atom1[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][1]; k++)
if (onethree[m][k] == angle_atom3[i][j]) {
flag[m][k] = 1;
break;
}
}
m = atom->map(angle_atom3[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][1]; k++)
if (onethree[m][k] == angle_atom1[i][j]) {
flag[m][k] = 1;
break;
}
}
}
}
if (num_dihedral) {
for (j = 0; j < num_dihedral[i]; j++) {
if (tag[i] != dihedral_atom2[i][j]) continue;
m = atom->map(dihedral_atom1[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][1]; k++)
if (onethree[m][k] == dihedral_atom3[i][j]) {
flag[m][k] = 1;
break;
}
}
m = atom->map(dihedral_atom3[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][1]; k++)
if (onethree[m][k] == dihedral_atom1[i][j]) {
flag[m][k] = 1;
break;
}
}
m = atom->map(dihedral_atom4[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][1]; k++)
if (onethree[m][k] == dihedral_atom2[i][j]) {
flag[m][k] = 1;
break;
}
}
}
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
for (k = 0; k < nspecial[i][1]; k++)
if (onethree[i][k] == outbuf[m].partnerID) {
flag[i][k] = 1;
break;
}
}
memory->destroy(outbuf);
for (i = 0; i < nlocal; i++) {
j = 0;
while (j < nspecial[i][1]) {
if (flag[i][j] == 0) {
onethree[i][j] = onethree[i][nspecial[i][1]-1];
flag[i][j] = flag[i][nspecial[i][1]-1];
nspecial[i][1]--;
} else j++;
}
}
memory->destroy(flag);
} else {
for (i = 0; i < nlocal; i++) nspecial[i][1] = 0;
}
onethreecount = 0.0;
for (i = 0; i < nlocal; i++) onethreecount += nspecial[i][1];
MPI_Allreduce(&onethreecount,&allcount,1,MPI_DOUBLE,MPI_SUM,world);
if (me == 0) {
if (screen)
fprintf(screen,
" %g = # of 1-3 neighbors after angle trim\n",allcount);
if (logfile)
fprintf(logfile,
" %g = # of 1-3 neighbors after angle trim\n",allcount);
}
}
void Special::dihedral_trim()
{
int i,j,k,m;
int *num_dihedral = atom->num_dihedral;
tagint **dihedral_atom1 = atom->dihedral_atom1;
tagint **dihedral_atom2 = atom->dihedral_atom2;
tagint **dihedral_atom4 = atom->dihedral_atom4;
int **nspecial = atom->nspecial;
tagint *tag = atom->tag;
int nlocal = atom->nlocal;
double onefourcount = 0.0;
for (i = 0; i < nlocal; i++) onefourcount += nspecial[i][2];
double allcount;
MPI_Allreduce(&onefourcount,&allcount,1,MPI_DOUBLE,MPI_SUM,world);
if (me == 0) {
if (screen)
fprintf(screen,
" %g = # of 1-4 neighbors before dihedral trim\n",allcount);
if (logfile)
fprintf(logfile,
" %g = # of 1-4 neighbors before dihedral trim\n",allcount);
}
if (num_dihedral && atom->ndihedrals) {
int nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_dihedral[i]; j++) {
if (tag[i] != dihedral_atom2[i][j]) continue;
m = atom->map(dihedral_atom1[i][j]);
if (m < 0 || m >= nlocal) nsend++;
m = atom->map(dihedral_atom4[i][j]);
if (m < 0 || m >= nlocal) nsend++;
}
}
int *proclist;
memory->create(proclist,nsend,"special:proclist");
PairRvous *inbuf = (PairRvous *)
memory->smalloc((bigint) nsend*sizeof(PairRvous),"special:inbuf");
nsend = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_dihedral[i]; j++) {
if (tag[i] != dihedral_atom2[i][j]) continue;
m = atom->map(dihedral_atom1[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = dihedral_atom1[i][j] % nprocs;
inbuf[nsend].atomID = dihedral_atom1[i][j];
inbuf[nsend].partnerID = dihedral_atom4[i][j];
nsend++;
}
m = atom->map(dihedral_atom4[i][j]);
if (m < 0 || m >= nlocal) {
proclist[nsend] = dihedral_atom4[i][j] % nprocs;
inbuf[nsend].atomID = dihedral_atom4[i][j];
inbuf[nsend].partnerID = dihedral_atom1[i][j];
nsend++;
}
}
}
char *buf;
int nreturn = comm->rendezvous(RVOUS,nsend,(char *) inbuf,sizeof(PairRvous),
0,proclist,
rendezvous_pairs,0,buf,sizeof(PairRvous),
(void *) this);
PairRvous *outbuf = (PairRvous *) buf;
memory->destroy(proclist);
memory->sfree(inbuf);
int max = 0;
for (i = 0; i < nlocal; i++)
max = MAX(max,nspecial[i][2]);
MPI_Allreduce(&max,&maxall,1,MPI_INT,MPI_MAX,world);
int **flag;
memory->create(flag,nlocal,maxall,"special:flag");
for (i = 0; i < nlocal; i++)
for (j = 0; j < nspecial[i][2]; j++)
flag[i][j] = 0;
for (i = 0; i < nlocal; i++) {
for (j = 0; j < num_dihedral[i]; j++) {
if (tag[i] != dihedral_atom2[i][j]) continue;
m = atom->map(dihedral_atom1[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][2]; k++)
if (onefour[m][k] == dihedral_atom4[i][j]) {
flag[m][k] = 1;
break;
}
}
m = atom->map(dihedral_atom4[i][j]);
if (m >= 0 && m < nlocal) {
for (k = 0; k < nspecial[m][2]; k++)
if (onefour[m][k] == dihedral_atom1[i][j]) {
flag[m][k] = 1;
break;
}
}
}
}
for (m = 0; m < nreturn; m++) {
i = atom->map(outbuf[m].atomID);
for (k = 0; k < nspecial[i][2]; k++)
if (onefour[i][k] == outbuf[m].partnerID) {
flag[i][k] = 1;
break;
}
}
memory->destroy(outbuf);
for (i = 0; i < nlocal; i++) {
j = 0;
while (j < nspecial[i][2]) {
if (flag[i][j] == 0) {
onefour[i][j] = onefour[i][nspecial[i][2]-1];
flag[i][j] = flag[i][nspecial[i][2]-1];
nspecial[i][2]--;
} else j++;
}
}
memory->destroy(flag);
} else {
for (i = 0; i < nlocal; i++) nspecial[i][2] = 0;
}
onefourcount = 0.0;
for (i = 0; i < nlocal; i++) onefourcount += nspecial[i][2];
MPI_Allreduce(&onefourcount,&allcount,1,MPI_DOUBLE,MPI_SUM,world);
if (me == 0) {
if (screen)
fprintf(screen,
" %g = # of 1-4 neighbors after dihedral trim\n",allcount);
if (logfile)
fprintf(logfile,
" %g = # of 1-4 neighbors after dihedral trim\n",allcount);
}
}
int Special::rendezvous_ids(int n, char *inbuf,
int &flag, int *& , char *& ,
void *ptr)
{
Special *sptr = (Special *) ptr;
Memory *memory = sptr->memory;
int *procowner;
tagint *atomIDs;
memory->create(procowner,n,"special:procowner");
memory->create(atomIDs,n,"special:atomIDs");
IDRvous *in = (IDRvous *) inbuf;
for (int i = 0; i < n; i++) {
procowner[i] = in[i].me;
atomIDs[i] = in[i].atomID;
}
sptr->nrvous = n;
sptr->procowner = procowner;
sptr->atomIDs = atomIDs;
flag = 0;
return 0;
}
int Special::rendezvous_pairs(int n, char *inbuf,
int &flag, int *&proclist, char *&outbuf,
void *ptr)
{
Special *sptr = (Special *) ptr;
Atom *atom = sptr->atom;
Memory *memory = sptr->memory;
atom->map_clear();
int nrvous = sptr->nrvous;
tagint *atomIDs = sptr->atomIDs;
for (int i = 0; i < nrvous; i++)
atom->map_one(atomIDs[i],i);
PairRvous *in = (PairRvous *) inbuf;
int *procowner = sptr->procowner;
memory->create(proclist,n,"special:proclist");
int m;
for (int i = 0; i < n; i++) {
m = atom->map(in[i].atomID);
proclist[i] = procowner[m];
}
outbuf = inbuf;
atom->map_init(0);
atom->nghost = 0;
atom->map_set();
flag = 1;
return n;
}
void Special::fix_alteration()
{
for (int ifix = 0; ifix < modify->nfix; ifix++)
if (modify->fix[ifix]->special_alter_flag)
modify->fix[ifix]->rebuild_special();
}
void Special::timer_output(double time1)
{
double time2 = MPI_Wtime();
if (comm->me == 0) {
if (screen) fprintf(screen," special bonds CPU = %g secs\n",time2-time1);
if (logfile) fprintf(logfile," special bonds CPU = %g secs\n",time2-time1);
}
}