/* ---------------------------------------------------------------------- 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 "compute_ke_molecule.h" #include "atom.h" #include "update.h" #include "force.h" #include "domain.h" #include "memory.h" #include "error.h" using namespace LAMMPS_NS; /* ---------------------------------------------------------------------- */ ComputeKEMolecule::ComputeKEMolecule(LAMMPS *lmp, int narg, char **arg) : Compute(lmp, narg, arg) { if (narg != 3) error->all("Illegal compute ke/molecule command"); if (atom->molecular == 0) error->all("Compute ke/molecule requires molecular atom style"); vector_flag = 1; extvector = 0; // setup molecule-based data nmolecules = molecules_in_group(idlo,idhi); size_vector = nmolecules; ke = (double *) memory->smalloc(nmolecules*sizeof(double), "ke/molecule:ke"); keall = (double *) memory->smalloc(nmolecules*sizeof(double), "ke/molecule:keall"); vector = keall; } /* ---------------------------------------------------------------------- */ ComputeKEMolecule::~ComputeKEMolecule() { memory->sfree(ke); memory->sfree(keall); } /* ---------------------------------------------------------------------- */ void ComputeKEMolecule::init() { int ntmp = molecules_in_group(idlo,idhi); if (ntmp != nmolecules) error->all("Molecule count changed in compute ke/molecule"); pfactor = 0.5 * force->mvv2e; } /* ---------------------------------------------------------------------- */ void ComputeKEMolecule::compute_vector() { int i,imol; double massone; invoked_array = update->ntimestep; int *mask = atom->mask; int *molecule = atom->molecule; int *type = atom->type; double **v = atom->v; double *rmass = atom->rmass; double *mass = atom->mass; int nlocal = atom->nlocal; for (i = 0; i < nmolecules; i++) ke[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--; if (rmass) massone = rmass[i]; else massone = mass[type[i]]; ke[imol] += massone * (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]); } MPI_Allreduce(ke,keall,nmolecules,MPI_DOUBLE,MPI_SUM,world); for (i = 0; i < nmolecules; i++) keall[i] = ke[i] * pfactor; }