/* ---------------------------------------------------------------------- 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. ------------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- Contributing author: Andres Jaramillo-Botero ------------------------------------------------------------------------- */ #include "mpi.h" #include "math.h" #include "compute_ke_eff.h" #include "atom.h" #include "update.h" #include "force.h" #include "domain.h" #include "group.h" #include "error.h" using namespace LAMMPS_NS; /* ---------------------------------------------------------------------- */ ComputeKEEff::ComputeKEEff(LAMMPS *lmp, int narg, char **arg) : Compute(lmp, narg, arg) { if (narg != 3) error->all("Illegal compute ke/eff command"); scalar_flag = 1; extscalar = 1; // error check if (!atom->ervel_flag) error->all("Compute ke/eff requires atom attribute ervel"); } /* ---------------------------------------------------------------------- */ void ComputeKEEff::init() { pfactor = 0.5 * force->mvv2e; } /* ---------------------------------------------------------------------- */ double ComputeKEEff::compute_scalar() { invoked_scalar = update->ntimestep; double **v = atom->v; double *ervel = atom->ervel; double *mass = atom->mass; int *spin = atom->spin; int *mask = atom->mask; int *type = atom->type; int nlocal = atom->nlocal; double ke = 0.0; if (mass) { for (int i = 0; i < nlocal; i++) if (mask[i] & groupbit) { ke += mass[type[i]]*(v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]); if (fabs(spin[i])==1) ke += 0.75*mass[type[i]]*ervel[i]*ervel[i]; } } MPI_Allreduce(&ke,&scalar,1,MPI_DOUBLE,MPI_SUM,world); scalar *= pfactor; return scalar; }