diff --git a/src/compute_displace_atom.cpp b/src/compute_displace_atom.cpp index ddeeb73aad..f17f54f8f8 100644 --- a/src/compute_displace_atom.cpp +++ b/src/compute_displace_atom.cpp @@ -71,7 +71,6 @@ ComputeDisplaceAtom::~ComputeDisplaceAtom() void ComputeDisplaceAtom::init() { // set fix which stores original atom coords - // check if is correct style int ifix = modify->find_fix(id_fix); if (ifix < 0) error->all("Could not find compute displace/atom fix ID"); @@ -109,6 +108,7 @@ void ComputeDisplaceAtom::compute_peratom() double xprd = domain->xprd; double yprd = domain->yprd; double zprd = domain->zprd; + int xbox,ybox,zbox; double dx,dy,dz; diff --git a/src/compute_msd.cpp b/src/compute_msd.cpp new file mode 100644 index 0000000000..764fcbd3f8 --- /dev/null +++ b/src/compute_msd.cpp @@ -0,0 +1,196 @@ +/* ---------------------------------------------------------------------- + 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 "string.h" +#include "compute_msd.h" +#include "atom.h" +#include "update.h" +#include "group.h" +#include "domain.h" +#include "modify.h" +#include "fix.h" +#include "error.h" + +using namespace LAMMPS_NS; + +/* ---------------------------------------------------------------------- */ + +ComputeMSD::ComputeMSD(LAMMPS *lmp, int narg, char **arg) : + Compute(lmp, narg, arg) +{ + if (narg < 3) error->all("Illegal compute msd command"); + + vector_flag = 1; + size_vector = 4; + extvector = 0; + + // optional args + + comflag = 0; + + int iarg = 3; + while (iarg < narg) { + if (strcmp(arg[iarg],"com") == 0) { + if (iarg+2 > narg) error->all("Illegal compute msd command"); + if (strcmp(arg[iarg+1],"no") == 0) comflag = 0; + else if (strcmp(arg[iarg+1],"yes") == 0) comflag = 1; + else error->all("Illegal compute msd command"); + iarg += 2; + } else error->all("Illegal compute msd command"); + } + + // create a new fix coord/original style with or without com keyword + // id = compute-ID + coord_original, fix group = compute group + + int n = strlen(id) + strlen("_coord_original") + 1; + id_fix = new char[n]; + strcpy(id_fix,id); + strcat(id_fix,"_coord_original"); + + char **newarg = new char*[5]; + newarg[0] = id_fix; + newarg[1] = group->names[igroup]; + newarg[2] = (char *) "coord/original"; + newarg[3] = (char *) "com"; + newarg[4] = (char *) "yes"; + if (comflag) modify->add_fix(5,newarg); + else modify->add_fix(3,newarg); + delete [] newarg; + + vector = new double[4]; +} + +/* ---------------------------------------------------------------------- */ + +ComputeMSD::~ComputeMSD() +{ + // check nfix in case all fixes have already been deleted + + if (modify->nfix) modify->delete_fix(id_fix); + + delete [] id_fix; + delete [] vector; +} + +/* ---------------------------------------------------------------------- */ + +void ComputeMSD::init() +{ + // set fix which stores original atom coords + + int ifix = modify->find_fix(id_fix); + if (ifix < 0) error->all("Could not find compute msd fix ID"); + fix = modify->fix[ifix]; + + // nmsd = # of atoms in group + + int *mask = atom->mask; + int nlocal = atom->nlocal; + + nmsd = 0; + for (int i = 0; i < nlocal; i++) + if (mask[i] & groupbit) nmsd++; + + int nmsd_all; + MPI_Allreduce(&nmsd,&nmsd_all,1,MPI_INT,MPI_SUM,world); + nmsd = nmsd_all; + + masstotal = group->mass(igroup); +} + +/* ---------------------------------------------------------------------- */ + +void ComputeMSD::compute_vector() +{ + invoked_vector = update->ntimestep; + + // cm = current center of mass + + double cm[3]; + if (comflag) group->xcm(igroup,masstotal,cm); + + // dx,dy,dz = displacement of atom from original position + // original unwrapped position is stored by fix + // relative to center of mass if comflag is set + // for triclinic, need to unwrap current atom coord via h matrix + + double **xoriginal = fix->vector_atom; + + double **x = atom->x; + int *mask = atom->mask; + int *image = atom->image; + int nlocal = atom->nlocal; + + double *h = domain->h; + double xprd = domain->xprd; + double yprd = domain->yprd; + double zprd = domain->zprd; + + double dx,dy,dz; + int xbox,ybox,zbox; + + double msd[4]; + msd[0] = msd[1] = msd[2] = msd[3] = 0.0; + + if (domain->triclinic == 0) { + for (int i = 0; i < nlocal; i++) + if (mask[i] & groupbit) { + xbox = (image[i] & 1023) - 512; + ybox = (image[i] >> 10 & 1023) - 512; + zbox = (image[i] >> 20) - 512; + if (comflag) { + dx = x[i][0] + xbox*xprd - cm[0] - xoriginal[i][0]; + dy = x[i][1] + ybox*yprd - cm[1] - xoriginal[i][1]; + dz = x[i][2] + zbox*zprd - cm[2] - xoriginal[i][2]; + } else { + dx = x[i][0] + xbox*xprd - xoriginal[i][0]; + dy = x[i][1] + ybox*yprd - xoriginal[i][1]; + dz = x[i][2] + zbox*zprd - xoriginal[i][2]; + } + msd[0] += dx*dx; + msd[1] += dy*dy; + msd[2] += dz*dz; + msd[3] += dx*dx + dy*dy + dz*dz; + } + + } else { + for (int i = 0; i < nlocal; i++) + if (mask[i] & groupbit) { + xbox = (image[i] & 1023) - 512; + ybox = (image[i] >> 10 & 1023) - 512; + zbox = (image[i] >> 20) - 512; + if (comflag) { + dx = x[i][0] + h[0]*xbox + h[5]*ybox + h[4]*zbox - + cm[0] - xoriginal[i][0]; + dy = x[i][1] + h[1]*ybox + h[3]*zbox - cm[1] - xoriginal[i][1]; + dz = x[i][2] + h[2]*zbox - cm[2] - xoriginal[i][2]; + } else { + dx = x[i][0] + h[0]*xbox + h[5]*ybox + h[4]*zbox - xoriginal[i][0]; + dy = x[i][1] + h[1]*ybox + h[3]*zbox - xoriginal[i][1]; + dz = x[i][2] + h[2]*zbox - xoriginal[i][2]; + } + msd[0] += dx*dx; + msd[1] += dy*dy; + msd[2] += dz*dz; + msd[3] += dx*dx + dy*dy + dz*dz; + } + } + + MPI_Allreduce(msd,vector,4,MPI_DOUBLE,MPI_SUM,world); + if (nmsd) { + vector[0] /= nmsd; + vector[1] /= nmsd; + vector[2] /= nmsd; + vector[3] /= nmsd; + } +} diff --git a/src/compute_msd.h b/src/compute_msd.h new file mode 100644 index 0000000000..37829e4aa2 --- /dev/null +++ b/src/compute_msd.h @@ -0,0 +1,37 @@ +/* ---------------------------------------------------------------------- + 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_MSD_H +#define COMPUTE_MSD_H + +#include "compute.h" + +namespace LAMMPS_NS { + +class ComputeMSD : public Compute { + public: + ComputeMSD(class LAMMPS *, int, char **); + ~ComputeMSD(); + void init(); + void compute_vector(); + + private: + int comflag,nmsd; + double masstotal; + char *id_fix; + class Fix *fix; +}; + +} + +#endif diff --git a/src/fix_coord_original.cpp b/src/fix_coord_original.cpp index c530d469c3..979c6a7c32 100644 --- a/src/fix_coord_original.cpp +++ b/src/fix_coord_original.cpp @@ -27,13 +27,28 @@ using namespace LAMMPS_NS; FixCoordOriginal::FixCoordOriginal(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg) { - if (narg != 3) error->all("Illegal fix coord/original command"); + if (narg < 3) error->all("Illegal fix coord/original command"); restart_peratom = 1; peratom_flag = 1; size_peratom = 3; peratom_freq = 1; + // optional args + + int comflag = 0; + + int iarg = 3; + while (iarg < narg) { + if (strcmp(arg[iarg],"com") == 0) { + if (iarg+2 > narg) error->all("Illegal fix coord/original command"); + if (strcmp(arg[iarg+1],"no") == 0) comflag = 0; + else if (strcmp(arg[iarg+1],"yes") == 0) comflag = 1; + else error->all("Illegal fix coord/original command"); + iarg += 2; + } else error->all("Illegal fix coord/original command"); + } + // perform initial allocation of atom-based array // register with Atom class @@ -42,16 +57,31 @@ FixCoordOriginal::FixCoordOriginal(LAMMPS *lmp, int narg, char **arg) : atom->add_callback(0); atom->add_callback(1); + // cm = original center of mass + + double cm[3]; + if (comflag) { + double masstotal = group->mass(igroup); + group->xcm(igroup,masstotal,cm); + } + // xoriginal = initial unwrapped positions of atoms - + // relative to center of mass if comflag is set + double **x = atom->x; int *mask = atom->mask; int *image = atom->image; int nlocal = atom->nlocal; for (int i = 0; i < nlocal; i++) { - if (mask[i] & groupbit) domain->unmap(x[i],image[i],xoriginal[i]); - else xoriginal[i][0] = xoriginal[i][1] = xoriginal[i][2] = 0.0; + if (mask[i] & groupbit) { + domain->unmap(x[i],image[i],xoriginal[i]); + if (comflag) { + xoriginal[i][0] -= cm[0]; + xoriginal[i][1] -= cm[1]; + xoriginal[i][2] -= cm[2]; + } + } else xoriginal[i][0] = xoriginal[i][1] = xoriginal[i][2] = 0.0; } } diff --git a/src/style.h b/src/style.h index 01a676fb72..fb7e3cf016 100644 --- a/src/style.h +++ b/src/style.h @@ -84,6 +84,7 @@ CommandStyle(write_restart,WriteRestart) #include "compute_heat_flux.h" #include "compute_ke.h" #include "compute_ke_atom.h" +#include "compute_msd.h" #include "compute_pe.h" #include "compute_pe_atom.h" #include "compute_pressure.h" @@ -110,6 +111,7 @@ ComputeStyle(group/group,ComputeGroupGroup) ComputeStyle(heat/flux,ComputeHeatFlux) ComputeStyle(ke,ComputeKE) ComputeStyle(ke/atom,ComputeKEAtom) +ComputeStyle(msd,ComputeMSD) ComputeStyle(pe,ComputePE) ComputeStyle(pe/atom,ComputePEAtom) ComputeStyle(pressure,ComputePressure)