Page MenuHomec4science

compute_pe_tally.cpp
No OneTemporary

File Metadata

Created
Mon, Nov 11, 10:42

compute_pe_tally.cpp

/* ----------------------------------------------------------------------
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 <string.h>
#include "compute_pe_tally.h"
#include "atom.h"
#include "group.h"
#include "pair.h"
#include "update.h"
#include "memory.h"
#include "error.h"
#include "force.h"
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
ComputePETally::ComputePETally(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg < 4) error->all(FLERR,"Illegal compute pe/tally command");
igroup2 = group->find(arg[3]);
if (igroup2 == -1)
error->all(FLERR,"Could not find compute pe/tally second group ID");
groupbit2 = group->bitmask[igroup2];
scalar_flag = 1;
vector_flag = 0;
peratom_flag = 1;
timeflag = 1;
comm_reverse = size_peratom_cols = 2;
extscalar = 1;
peflag = 1; // we need Pair::ev_tally() to be run
did_compute = invoked_peratom = invoked_scalar = -1;
nmax = -1;
eatom = NULL;
vector = new double[size_peratom_cols];
}
/* ---------------------------------------------------------------------- */
ComputePETally::~ComputePETally()
{
if (force && force->pair) force->pair->del_tally_callback(this);
memory->destroy(eatom);
delete[] vector;
}
/* ---------------------------------------------------------------------- */
void ComputePETally::init()
{
if (force->pair == NULL)
error->all(FLERR,"Trying to use compute pe/tally with no pair style");
else
force->pair->add_tally_callback(this);
if (force->pair->single_enable == 0 || force->pair->manybody_flag)
error->all(FLERR,"Compute pe/tally used with incompatible pair style.");
if ((comm->me == 0) && (force->bond || force->angle || force->dihedral
|| force->improper || force->kspace))
error->warning(FLERR,"Compute pe/tally only called from pair style");
did_compute = -1;
}
/* ---------------------------------------------------------------------- */
void ComputePETally::pair_tally_callback(int i, int j, int nlocal, int newton,
double evdwl, double ecoul, double,
double, double, double)
{
const int ntotal = atom->nlocal + atom->nghost;
const int * const mask = atom->mask;
// do setup work that needs to be done only once per timestep
if (did_compute != update->ntimestep) {
did_compute = update->ntimestep;
// grow local eatom array if necessary
// needs to be atom->nmax in length
if (atom->nmax > nmax) {
memory->destroy(eatom);
nmax = atom->nmax;
memory->create(eatom,nmax,size_peratom_cols,"pe/tally:eatom");
array_atom = eatom;
}
// clear storage as needed
if (newton) {
for (int i=0; i < ntotal; ++i)
eatom[i][0] = eatom[i][1] = 0.0;
} else {
for (int i=0; i < atom->nlocal; ++i)
eatom[i][0] = eatom[i][1] = 0.0;
}
vector[0] = etotal[0] = vector[1] = etotal[1] = 0.0;
}
if ( ((mask[i] & groupbit) && (mask[j] & groupbit2))
|| ((mask[i] & groupbit2) && (mask[j] & groupbit)) ){
evdwl *= 0.5; ecoul *= 0.5;
if (newton || i < nlocal) {
etotal[0] += evdwl; eatom[i][0] += evdwl;
etotal[1] += ecoul; eatom[i][1] += ecoul;
}
if (newton || j < nlocal) {
etotal[0] += evdwl; eatom[j][0] += evdwl;
etotal[1] += ecoul; eatom[j][1] += ecoul;
}
}
}
/* ---------------------------------------------------------------------- */
int ComputePETally::pack_reverse_comm(int n, int first, double *buf)
{
int i,m,last;
m = 0;
last = first + n;
for (i = first; i < last; i++) {
buf[m++] = eatom[i][0];
buf[m++] = eatom[i][1];
}
return m;
}
/* ---------------------------------------------------------------------- */
void ComputePETally::unpack_reverse_comm(int n, int *list, double *buf)
{
int i,j,m;
m = 0;
for (i = 0; i < n; i++) {
j = list[i];
eatom[j][0] += buf[m++];
eatom[j][1] += buf[m++];
}
}
/* ---------------------------------------------------------------------- */
double ComputePETally::compute_scalar()
{
invoked_scalar = update->ntimestep;
if ((did_compute != invoked_scalar) || (update->eflag_global != invoked_scalar))
error->all(FLERR,"Energy was not tallied on needed timestep");
// sum accumulated energies across procs
MPI_Allreduce(etotal,vector,size_peratom_cols,MPI_DOUBLE,MPI_SUM,world);
scalar = vector[0]+vector[1];
return scalar;
}
/* ---------------------------------------------------------------------- */
void ComputePETally::compute_peratom()
{
invoked_peratom = update->ntimestep;
if ((did_compute != invoked_peratom) || (update->eflag_global != invoked_peratom))
error->all(FLERR,"Energy was not tallied on needed timestep");
// collect contributions from ghost atoms
if (force->newton_pair) {
comm->reverse_comm_compute(this);
const int nall = atom->nlocal + atom->nghost;
for (int i = atom->nlocal; i < nall; ++i)
eatom[i][0] = eatom[i][1] = 0.0;
}
}
/* ----------------------------------------------------------------------
memory usage of local atom-based array
------------------------------------------------------------------------- */
double ComputePETally::memory_usage()
{
double bytes = nmax*size_peratom_cols * sizeof(double);
return bytes;
}

Event Timeline