git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@1507 f3b2605a-c512-4ea7-a41b-209d697bcdaa

This commit is contained in:
sjplimp 2008-02-14 23:26:03 +00:00
parent 944bc255c3
commit ab1a5ae04f
3 changed files with 244 additions and 0 deletions

201
src/compute_group_group.cpp Normal file
View File

@ -0,0 +1,201 @@
/* ----------------------------------------------------------------------
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.
------------------------------------------------------------------------- */
#include "mpi.h"
#include "compute_group_group.h"
#include "atom.h"
#include "force.h"
#include "pair.h"
#include "neighbor.h"
#include "neigh_request.h"
#include "neigh_list.h"
#include "group.h"
#include "error.h"
using namespace LAMMPS_NS;
#define INVOKED_SCALAR 1
#define INVOKED_VECTOR 2
/* ---------------------------------------------------------------------- */
ComputeGroupGroup::ComputeGroupGroup(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg != 4) error->all("Illegal compute group/group command");
scalar_flag = vector_flag = 1;
size_vector = 3;
extscalar = 1;
extvector = 1;
jgroup = group->find(arg[3]);
if (jgroup == -1) error->all("Compute group/group group ID does not exist");
jgroupbit = group->bitmask[jgroup];
vector = new double[3];
}
/* ---------------------------------------------------------------------- */
ComputeGroupGroup::~ComputeGroupGroup()
{
delete [] vector;
}
/* ---------------------------------------------------------------------- */
void ComputeGroupGroup::init()
{
if (force->pair == NULL) error->all("Bad");
pair = force->pair;
cutsq = force->pair->cutsq;
// need an occasional half neighbor list
int irequest = neighbor->request((void *) this);
neighbor->requests[irequest]->pair = 0;
neighbor->requests[irequest]->compute = 1;
neighbor->requests[irequest]->occasional = 1;
}
/* ---------------------------------------------------------------------- */
void ComputeGroupGroup::init_list(int id, NeighList *ptr)
{
list = ptr;
}
/* ---------------------------------------------------------------------- */
double ComputeGroupGroup::compute_scalar()
{
invoked |= INVOKED_SCALAR;
invoked |= INVOKED_VECTOR;
interact();
return scalar;
}
/* ---------------------------------------------------------------------- */
void ComputeGroupGroup::compute_vector()
{
invoked |= INVOKED_SCALAR;
invoked |= INVOKED_VECTOR;
interact();
}
/* ---------------------------------------------------------------------- */
void ComputeGroupGroup::interact()
{
int i,j,ii,jj,inum,jnum,itype,jtype;
double xtmp,ytmp,ztmp,delx,dely,delz;
double rsq,eng,fpair,factor_coul,factor_lj;
int *ilist,*jlist,*numneigh,**firstneigh;
double **x = atom->x;
int *type = atom->type;
int *mask = atom->mask;
int nlocal = atom->nlocal;
int nall = nlocal + atom->nghost;
double *special_coul = force->special_coul;
double *special_lj = force->special_lj;
int newton_pair = force->newton_pair;
// invoke half 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;
// loop over neighbors of my atoms
// skip if i,j are not in 2 groups
double one[4],all[4];
one[0] = one[1] = one[2] = one[3] = 0.0;
for (ii = 0; ii < inum; ii++) {
i = ilist[ii];
if (mask[i] & groupbit) othergroupbit = jgroupbit;
else if (mask[i] & jgroupbit) othergroupbit = groupbit;
else continue;
xtmp = x[i][0];
ytmp = x[i][1];
ztmp = x[i][2];
itype = type[i];
jlist = firstneigh[i];
jnum = numneigh[i];
for (jj = 0; jj < jnum; jj++) {
j = jlist[jj];
if (j < nall) factor_coul = factor_lj = 1.0;
else {
factor_coul = special_coul[j/nall];
factor_lj = special_lj[j/nall];
j %= nall;
}
if (!(mask[j] & othergroupbit)) continue;
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
delz = ztmp - x[j][2];
rsq = delx*delx + dely*dely + delz*delz;
jtype = type[j];
if (rsq < cutsq[itype][jtype]) {
eng = pair->single(i,j,itype,jtype,rsq,factor_coul,factor_lj,fpair);
// energy only computed once so tally full amount
// force tally is jgroup acting on igroup
if (newton_pair || j < nlocal) {
one[0] += eng;
if (othergroupbit == jgroupbit) {
one[1] += delx*fpair;
one[2] += dely*fpair;
one[3] += delz*fpair;
}
if (othergroupbit == groupbit) {
one[1] -= delx*fpair;
one[2] -= dely*fpair;
one[3] -= delz*fpair;
}
// energy computed twice so tally half amount
// only tally force if I own igroup atom
} else {
one[0] += 0.5*eng;
if (othergroupbit == jgroupbit) {
one[1] += delx*fpair;
one[2] += dely*fpair;
one[3] += delz*fpair;
}
}
}
}
}
MPI_Allreduce(one,all,4,MPI_DOUBLE,MPI_SUM,world);
scalar = all[0];
vector[0] = all[1]; vector[1] = all[2]; vector[2] = all[3];
}

41
src/compute_group_group.h Normal file
View File

@ -0,0 +1,41 @@
/* ----------------------------------------------------------------------
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.
------------------------------------------------------------------------- */
#ifndef COMPUTE_GROUP_GROUP_H
#define COMPUTE_GROUP_GROUP_H
#include "compute.h"
namespace LAMMPS_NS {
class ComputeGroupGroup : public Compute {
public:
ComputeGroupGroup(class LAMMPS *, int, char **);
~ComputeGroupGroup();
void init();
void init_list(int, class NeighList *);
double compute_scalar();
void compute_vector();
private:
int jgroup,jgroupbit,othergroupbit;
double **cutsq;
class Pair *pair;
class NeighList *list;
void interact();
};
}
#endif

View File

@ -79,6 +79,7 @@ CommandStyle(write_restart,WriteRestart)
#include "compute_centro_atom.h"
#include "compute_coord_atom.h"
#include "compute_displace_atom.h"
#include "compute_group_group.h"
#include "compute_ke_atom.h"
#include "compute_pe.h"
#include "compute_pe_atom.h"
@ -98,6 +99,7 @@ CommandStyle(write_restart,WriteRestart)
ComputeStyle(centro/atom,ComputeCentroAtom)
ComputeStyle(coord/atom,ComputeCoordAtom)
ComputeStyle(displace/atom,ComputeDisplaceAtom)
ComputeStyle(group/group,ComputeGroupGroup)
ComputeStyle(ke/atom,ComputeKEAtom)
ComputeStyle(pe,ComputePE)
ComputeStyle(pe/atom,ComputePEAtom)