Commit 972f8a08 authored by sjplimp's avatar sjplimp
Browse files

git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@1639 f3b2605a-c512-4ea7-a41b-209d697bcdaa
parent fb1abedd
Loading
Loading
Loading
Loading

src/compute_ke.cpp

0 → 100644
+74 −0
Original line number Diff line number Diff line
/* ----------------------------------------------------------------------
   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 "mpi.h"
#include "compute_ke.h"
#include "atom.h"
#include "force.h"
#include "domain.h"
#include "group.h"
#include "error.h"

using namespace LAMMPS_NS;

#define INVOKED_SCALAR 1

/* ---------------------------------------------------------------------- */

ComputeKE::ComputeKE(LAMMPS *lmp, int narg, char **arg) :
  Compute(lmp, narg, arg)
{
  if (narg != 3) error->all("Illegal compute ke command");

  scalar_flag = 1;
  extscalar = 1;
}

/* ---------------------------------------------------------------------- */

void ComputeKE::init()
{
  pfactor = 0.5 * force->mvv2e;
}

/* ---------------------------------------------------------------------- */

double ComputeKE::compute_scalar()
{
  invoked |= INVOKED_SCALAR;

  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 (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]);
	  
  } else {
    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]);
  }

  MPI_Allreduce(&ke,&scalar,1,MPI_DOUBLE,MPI_SUM,world);
  scalar *= pfactor;
  return scalar;
}

src/compute_ke.h

0 → 100644
+33 −0
Original line number Diff line number Diff line
/* ----------------------------------------------------------------------
   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.
------------------------------------------------------------------------- */

#ifndef COMPUTE_KE_H
#define COMPUTE_KE_H

#include "compute.h"

namespace LAMMPS_NS {

class ComputeKE : public Compute {
 public:
  ComputeKE(class LAMMPS *, int, char **);
  void init();
  double compute_scalar();

 private:
  double pfactor;
};

}

#endif
+2 −0
Original line number Diff line number Diff line
@@ -80,6 +80,7 @@ CommandStyle(write_restart,WriteRestart)
#include "compute_coord_atom.h"
#include "compute_displace_atom.h"
#include "compute_group_group.h"
#include "compute_ke.h"
#include "compute_ke_atom.h"
#include "compute_pe.h"
#include "compute_pe_atom.h"
@@ -101,6 +102,7 @@ ComputeStyle(centro/atom,ComputeCentroAtom)
ComputeStyle(coord/atom,ComputeCoordAtom)
ComputeStyle(displace/atom,ComputeDisplaceAtom)
ComputeStyle(group/group,ComputeGroupGroup)
ComputeStyle(ke,ComputeKE)
ComputeStyle(ke/atom,ComputeKEAtom)
ComputeStyle(pe,ComputePE)
ComputeStyle(pe/atom,ComputePEAtom)