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
/* ----------------------------------------------------------------------
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 "compute_ke.h"
#include <mpi.h>
#include "atom.h"
#include "update.h"
#include "force.h"
#include "error.h"
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
ComputeKE::ComputeKE(LAMMPS *lmp, int narg, char **arg) :
Compute(lmp, narg, arg)
{
if (narg != 3) error->all(FLERR,"Illegal compute ke command");
scalar_flag = 1;
extscalar = 1;
}
/* ---------------------------------------------------------------------- */
void ComputeKE::init()
{
pfactor = 0.5 * force->mvv2e;
}
/* ---------------------------------------------------------------------- */
double ComputeKE::compute_scalar()
{
invoked_scalar = update->ntimestep;
double **v = atom->v;
double *rmass = atom->rmass;
double *mass = atom->mass;
int *mask = atom->mask;
int *type = atom->type;
int nlocal = atom->nlocal;
double ke = 0.0;
if (rmass) {
for (int i = 0; i < nlocal; i++)
if (mask[i] & groupbit)
ke += rmass[i] * (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]);
} else {
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]);
}
MPI_Allreduce(&ke,&scalar,1,MPI_DOUBLE,MPI_SUM,world);
scalar *= pfactor;
return scalar;
}