/* ---------------------------------------------------------------------- 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 "fix_my_momrot_chunk.h" #include #include #include "atom.h" #include "domain.h" #include "update.h" #include "force.h" #include "modify.h" #include "compute_chunk_atom.h" #include "compute_com_chunk.h" #include "compute_vcm_chunk.h" #include "compute_omega_chunk.h" #include "memory.h" #include "error.h" using namespace LAMMPS_NS; using namespace FixConst; #define SMALL 1.0e-10 /* ---------------------------------------------------------------------- */ FixMomentumChunk::FixMomentumChunk(LAMMPS *lmp, int narg, char **arg) : Fix(lmp, narg, arg), idchunk(NULL), idcom(NULL),idvcm(NULL),idomega(NULL) { if (narg < 9) error->all(FLERR,"Illegal fix momentum/chun command"); nevery = force->inumeric(FLERR,arg[3]); if (nevery <= 0) error->all(FLERR,"Illegal fix momentum/chun command"); linear = angular= 0; int n = strlen(arg[4]) + 1; idchunk = new char[n]; strcpy(idchunk,arg[4]); n = strlen(arg[5]) + 1; idcom = new char[n]; strcpy(idcom,arg[5]); n = strlen(arg[6]) + 1; idvcm = new char[n]; strcpy(idvcm,arg[6]); n = strlen(arg[7]) + 1; idomega = new char[n]; strcpy(idomega,arg[7]); int iarg = 8; while (iarg < narg) { if (strcmp(arg[iarg],"linear") == 0) { linear = 1; iarg += 1; } else if (strcmp(arg[iarg],"angular") == 0) { angular = 1; iarg += 1; } else error->all(FLERR,"Illegal fix momentum/chunk command"); } if (linear == 0 && angular == 0) error->all(FLERR,"Illegal fix momentum/chunk command"); nchunk = 0; } /* ---------------------------------------------------------------------- */ FixMomentumChunk::~FixMomentumChunk() { // decrement lock counter in compute chunk/atom, it if still exists int icompute = modify->find_compute(idchunk); if (icompute >= 0) { cchunk = (ComputeChunkAtom *) modify->compute[icompute]; cchunk->unlock(this); cchunk->lockcount--; } delete [] idchunk; delete [] idcom; } /* ---------------------------------------------------------------------- */ int FixMomentumChunk::setmask() { int mask = 0; mask |= END_OF_STEP; return mask; } /* ---------------------------------------------------------------------- */ void FixMomentumChunk::init() { // current indices for idchunk and idcom int icompute = modify->find_compute(idchunk); if (icompute < 0) error->all(FLERR,"Chunk/atom compute does not exist for fix momentum/chunk"); cchunk = (ComputeChunkAtom *) modify->compute[icompute]; if (strcmp(cchunk->style,"chunk/atom") != 0) error->all(FLERR,"Fix momentum/chunk does not use chunk/atom compute"); icompute = modify->find_compute(idcom); if (icompute < 0) error->all(FLERR,"Com/chunk compute does not exist for fix momentum/chunk"); ccom = (ComputeCOMChunk *) modify->compute[icompute]; if (strcmp(ccom->style,"com/chunk") != 0) error->all(FLERR,"Fix momentum/chunk does not use com/chunk compute"); icompute = modify->find_compute(idvcm); if (icompute < 0) error->all(FLERR,"vcm/chunk compute does not exist for fix momentum/chunk"); cvcm = (ComputeVCMChunk *) modify->compute[icompute]; if (strcmp(cvcm->style,"vcm/chunk") != 0) error->all(FLERR,"Fix momentum/chunk does not use vcm/chunk compute"); icompute = modify->find_compute(idomega); if (icompute < 0) error->all(FLERR,"omega/chunk compute does not exist for fix momentum/chunk"); comega = (ComputeOmegaChunk *) modify->compute[icompute]; if (strcmp(comega->style,"omega/chunk") != 0) error->all(FLERR,"Fix momentum/chunk does not use omega/chunk compute"); // check that idchunk is consistent with ccom->idchunk if (strcmp(idchunk,ccom->idchunk) != 0) error->all(FLERR,"Fix spring chunk chunkID not same as comID chunkID"); // check that idchunk is consistent with cvcm->idchunk //idchunk is private data of VCM ad OMEGA, no check for them // if (strcmp(idchunk,cvcm->idchunk) != 0) // error->all(FLERR,"Fix momentum chunk chunkID not same as vcmID chunkID"); // // check that idchunk is consistent with comega->idchunk // if (strcmp(idchunk,comega->idchunk) != 0) // error->all(FLERR,"Fix momentum chunk chunkID not same as omegaID chunkID"); } /* ---------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- */ void FixMomentumChunk::end_of_step() { int i,m; double dx,dy,dz; // calculate current centers of mass for each chunk // extract pointers from idchunk and idcom ccom->compute_array(); cvcm->compute_array(); comega->compute_array(); nchunk = cchunk->nchunk; int *ichunk = cchunk->ichunk; double *masstotal = ccom->masstotal; double **com = ccom->array; double **vcm = cvcm->array; double **omega = comega->array; // apply removing translational and rotational velocity to atoms in each chunk double **v = atom->v; int *mask = atom->mask; const int nlocal = atom->nlocal; int *type = atom->type; double *mass = atom->mass; double *rmass = atom->rmass; if (linear) { for (i = 0; i < nlocal; i++) { m = ichunk[i]-1; if (m < 0) continue; v[i][0] -= vcm[m][0]; v[i][1] -= vcm[m][1]; v[i][2] -= vcm[m][2]; } } if (angular) { // must use unwrapped coords to compute omega to v correctly double **x = atom->x; imageint *image = atom->image; double unwrap[3]; for (i = 0; i < nlocal; i++) { m = ichunk[i]-1; if (m < 0) continue; domain->unmap(x[i],image[i],unwrap); dx = unwrap[0] - com[m][0]; dy = unwrap[1] - com[m][1]; dz = unwrap[2] - com[m][2]; v[i][0] -= omega[m][1]*dz - omega[m][2]*dy; v[i][1] -= omega[m][2]*dx - omega[m][0]*dz; v[i][2] -= omega[m][0]*dy - omega[m][1]*dx; } } }