/* ---------------------------------------------------------------------- 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 "math.h" #include "string.h" #include "compute_pe_molecule.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 "domain.h" #include "memory.h" #include "error.h" using namespace LAMMPS_NS; /* ---------------------------------------------------------------------- */ ComputePEMolecule::ComputePEMolecule(LAMMPS *lmp, int narg, char **arg) : Compute(lmp, narg, arg) { if (narg < 3) error->all("Illegal compute pe/molecule command"); vector_flag = 1; extvector = 0; peratom_flag = 1; size_peratom_cols = 0; peatomflag = 1; timeflag = 1; comm_reverse = 1; nmolecules = molecules_in_group(idlo,idhi); size_vector = nmolecules; // setup molecule-based data pe = (double *) memory->smalloc(nmolecules*sizeof(double),"pe/molecule:pe"); peall = (double *) memory->smalloc(nmolecules*sizeof(double),"pe/molecule:peall"); vector = peall; // flags for potential types if (narg == 3) { pairflag = 1; bondflag = angleflag = dihedralflag = improperflag = 1; } else { pairflag = 0; bondflag = angleflag = dihedralflag = improperflag = 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 error->all("Illegal compute pe/molecule command"); iarg++; } } nmax = 0; energy = NULL; } /* ---------------------------------------------------------------------- */ ComputePEMolecule::~ComputePEMolecule() { memory->sfree(energy); memory->sfree(pe); memory->sfree(peall); } /* ---------------------------------------------------------------------- */ void ComputePEMolecule::compute_vector() { int i,imol; int *mask = atom->mask; int *molecule = atom->molecule; invoked_array = update->ntimestep; // grow local energy array if necessary // needs to be atom->nmax in length if (atom->nmax > nmax) { memory->sfree(energy); nmax = atom->nmax; energy = (double *) memory->smalloc(nmax*sizeof(double),"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); // zero energy of atoms not in group // only do this after comm since ghost contributions must be included for (i = 0; i < nlocal; i++) if (!(mask[i] & groupbit)) energy[i] = 0.0; for (i = 0; i < nmolecules; i++) pe[i] = 0.0; for (int i = 0; i < nlocal; i++) if (mask[i] & groupbit) { imol = molecule[i]; if (molmap) imol = molmap[imol-idlo]; else imol--; pe[imol] += energy[i]; } MPI_Allreduce(pe,peall,nmolecules,MPI_DOUBLE,MPI_SUM,world); } /* ---------------------------------------------------------------------- */ int ComputePEMolecule::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 ComputePEMolecule::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 ComputePEMolecule::memory_usage() { double bytes = nmax * sizeof(double); bytes += 2 * nmolecules * sizeof(double); if (molmap) bytes += (idhi-idlo+1); return bytes; }