1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193
|
/* ----------------------------------------------------------------------
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_atom.h"
#include "atom.h"
#include "update.h"
#include "comm.h"
#include "force.h"
#include "pair.h"
#include "bond.h"
#include "angle.h"
#include "dihedral.h"
#include "improper.h"
#include "kspace.h"
#include "memory.h"
#include "error.h"
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
ComputePEAtom::ComputePEAtom(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg < 3) error->all(FLERR,"Illegal compute pe/atom command");
peratom_flag = 1;
size_peratom_cols = 0;
peatomflag = 1;
timeflag = 1;
comm_reverse = 1;
if (narg == 3) {
pairflag = 1;
bondflag = angleflag = dihedralflag = improperflag = 1;
kspaceflag = 1;
} else {
pairflag = 0;
bondflag = angleflag = dihedralflag = improperflag = 0;
kspaceflag = 0;
int iarg = 3;
while (iarg < narg) {
if (strcmp(arg[iarg],"pair") == 0) pairflag = 1;
else if (strcmp(arg[iarg],"bond") == 0) bondflag = 1;
else if (strcmp(arg[iarg],"angle") == 0) angleflag = 1;
else if (strcmp(arg[iarg],"dihedral") == 0) dihedralflag = 1;
else if (strcmp(arg[iarg],"improper") == 0) improperflag = 1;
else if (strcmp(arg[iarg],"kspace") == 0) kspaceflag = 1;
else error->all(FLERR,"Illegal compute pe/atom command");
iarg++;
}
}
nmax = 0;
energy = NULL;
}
/* ---------------------------------------------------------------------- */
ComputePEAtom::~ComputePEAtom()
{
memory->destroy(energy);
}
/* ---------------------------------------------------------------------- */
void ComputePEAtom::compute_peratom()
{
int i;
invoked_peratom = update->ntimestep;
if (update->eflag_atom != invoked_peratom)
error->all(FLERR,"Per-atom energy was not tallied on needed timestep");
// grow local energy array if necessary
// needs to be atom->nmax in length
if (atom->nmax > nmax) {
memory->destroy(energy);
nmax = atom->nmax;
memory->create(energy,nmax,"pe/atom:energy");
vector_atom = energy;
}
// npair includes ghosts if either newton flag is set
// b/c some bonds/dihedrals call pair::ev_tally with pairwise info
// nbond includes ghosts if newton_bond is set
// ntotal includes ghosts if either newton flag is set
int nlocal = atom->nlocal;
int npair = nlocal;
int nbond = nlocal;
int ntotal = nlocal;
if (force->newton) npair += atom->nghost;
if (force->newton_bond) nbond += atom->nghost;
if (force->newton) ntotal += atom->nghost;
// clear local energy array
for (i = 0; i < ntotal; i++) energy[i] = 0.0;
// add in per-atom contributions from each force
if (pairflag && force->pair) {
double *eatom = force->pair->eatom;
for (i = 0; i < npair; i++) energy[i] += eatom[i];
}
if (bondflag && force->bond) {
double *eatom = force->bond->eatom;
for (i = 0; i < nbond; i++) energy[i] += eatom[i];
}
if (angleflag && force->angle) {
double *eatom = force->angle->eatom;
for (i = 0; i < nbond; i++) energy[i] += eatom[i];
}
if (dihedralflag && force->dihedral) {
double *eatom = force->dihedral->eatom;
for (i = 0; i < nbond; i++) energy[i] += eatom[i];
}
if (improperflag && force->improper) {
double *eatom = force->improper->eatom;
for (i = 0; i < nbond; i++) energy[i] += eatom[i];
}
// communicate ghost energy between neighbor procs
if (force->newton) comm->reverse_comm_compute(this);
// KSpace contribution is already per local atom
if (kspaceflag && force->kspace) {
double *eatom = force->kspace->eatom;
for (i = 0; i < nlocal; i++) energy[i] += eatom[i];
}
// zero energy of atoms not in group
// only do this after comm since ghost contributions must be included
int *mask = atom->mask;
for (i = 0; i < nlocal; i++)
if (!(mask[i] & groupbit)) energy[i] = 0.0;
}
/* ---------------------------------------------------------------------- */
int ComputePEAtom::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++] = energy[i];
return 1;
}
/* ---------------------------------------------------------------------- */
void ComputePEAtom::unpack_reverse_comm(int n, int *list, double *buf)
{
int i,j,m;
m = 0;
for (i = 0; i < n; i++) {
j = list[i];
energy[j] += buf[m++];
}
}
/* ----------------------------------------------------------------------
memory usage of local atom-based array
------------------------------------------------------------------------- */
double ComputePEAtom::memory_usage()
{
double bytes = nmax * sizeof(double);
return bytes;
}
|