/* ---------------------------------------------------------------------- This is the ██╗ ██╗ ██████╗ ██████╗ ██████╗ ██╗ ██╗████████╗███████╗ ██║ ██║██╔════╝ ██╔════╝ ██╔════╝ ██║ ██║╚══██╔══╝██╔════╝ ██║ ██║██║ ███╗██║ ███╗██║ ███╗███████║ ██║ ███████╗ ██║ ██║██║ ██║██║ ██║██║ ██║██╔══██║ ██║ ╚════██║ ███████╗██║╚██████╔╝╚██████╔╝╚██████╔╝██║ ██║ ██║ ███████║ ╚══════╝╚═╝ ╚═════╝ ╚═════╝ ╚═════╝ ╚═╝ ╚═╝ ╚═╝ ╚══════╝® DEM simulation engine, released by DCS Computing Gmbh, Linz, Austria http://www.dcs-computing.com, office@dcs-computing.com LIGGGHTS® is part of CFDEM®project: http://www.liggghts.com | http://www.cfdem.com Core developer and main author: Christoph Kloss, christoph.kloss@dcs-computing.com LIGGGHTS® is open-source, distributed under the terms of the GNU Public License, version 2 or later. It is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. You should have received a copy of the GNU General Public License along with LIGGGHTS®. If not, see http://www.gnu.org/licenses . See also top-level README and LICENSE files. LIGGGHTS® and CFDEM® are registered trade marks of DCS Computing GmbH, the producer of the LIGGGHTS® software and the CFDEM®coupling software See http://www.cfdem.com/terms-trademark-policy for details. ------------------------------------------------------------------------- Contributing author and copyright for this file: This file is from LAMMPS 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. ------------------------------------------------------------------------- */ #include #include "compute_erotate_asphere.h" #include "math_extra.h" #include "atom.h" #include "atom_vec_ellipsoid.h" #include "atom_vec_line.h" #include "atom_vec_tri.h" #include "update.h" #include "force.h" #include "memory.h" #include "error.h" using namespace LAMMPS_NS; /* ---------------------------------------------------------------------- */ ComputeERotateAsphere::ComputeERotateAsphere(LAMMPS *lmp, int & iarg, int narg, char **arg) : Compute(lmp, iarg, narg, arg) { if (narg != iarg) error->all(FLERR,"Illegal compute erotate/asphere command"); scalar_flag = 1; extscalar = 1; } /* ---------------------------------------------------------------------- */ void ComputeERotateAsphere::init() { // error check avec_ellipsoid = (AtomVecEllipsoid *) atom->style_match("ellipsoid"); avec_line = (AtomVecLine *) atom->style_match("line"); avec_tri = (AtomVecTri *) atom->style_match("tri"); if (!avec_ellipsoid && !avec_line && !avec_tri) error->all(FLERR,"Compute erotate/asphere requires " "atom style ellipsoid or line or tri"); // check that all particles are finite-size // no point particles allowed, spherical is OK int *ellipsoid = atom->ellipsoid; int *line = atom->line; int *tri = atom->tri; int *mask = atom->mask; int nlocal = atom->nlocal; int flag; for (int i = 0; i < nlocal; i++) if (mask[i] & groupbit) { flag = 0; if (ellipsoid && ellipsoid[i] >= 0) flag = 1; if (line && line[i] >= 0) flag = 1; if (tri && tri[i] >= 0) flag = 1; if (!flag) error->one(FLERR,"Compute erotate/asphere requires extended particles"); } pfactor = 0.5 * force->mvv2e; } /* ---------------------------------------------------------------------- */ double ComputeERotateAsphere::compute_scalar() { invoked_scalar = update->ntimestep; AtomVecEllipsoid::Bonus *ebonus = NULL; if (avec_ellipsoid) ebonus = avec_ellipsoid->bonus; AtomVecLine::Bonus *lbonus = NULL; if (avec_line) lbonus = avec_line->bonus; AtomVecTri::Bonus *tbonus = NULL; if (avec_tri) tbonus = avec_tri->bonus; int *ellipsoid = atom->ellipsoid; int *line = atom->line; int *tri = atom->tri; double **omega = atom->omega; double **angmom = atom->angmom; double *rmass = atom->rmass; int *mask = atom->mask; int nlocal = atom->nlocal; // sum rotational energy for each particle // no point particles since divide by inertia double length; double *shape,*quat; double wbody[3],inertia[3]; double rot[3][3]; double erotate = 0.0; for (int i = 0; i < nlocal; i++) if (mask[i] & groupbit) { if (ellipsoid && ellipsoid[i] >= 0) { shape = ebonus[ellipsoid[i]].shape; quat = ebonus[ellipsoid[i]].quat; // principal moments of inertia inertia[0] = rmass[i] * (shape[1]*shape[1]+shape[2]*shape[2]) / 5.0; inertia[1] = rmass[i] * (shape[0]*shape[0]+shape[2]*shape[2]) / 5.0; inertia[2] = rmass[i] * (shape[0]*shape[0]+shape[1]*shape[1]) / 5.0; // wbody = angular velocity in body frame MathExtra::quat_to_mat(quat,rot); MathExtra::transpose_matvec(rot,angmom[i],wbody); wbody[0] /= inertia[0]; wbody[1] /= inertia[1]; wbody[2] /= inertia[2]; erotate += inertia[0]*wbody[0]*wbody[0] + inertia[1]*wbody[1]*wbody[1] + inertia[2]*wbody[2]*wbody[2]; } else if (line && line[i] >= 0) { length = lbonus[line[i]].length; erotate += (omega[i][0]*omega[i][0] + omega[i][1]*omega[i][1] + omega[i][2]*omega[i][2]) * length*length*rmass[i] / 12.0; } else if (tri && tri[i] >= 0) { // principal moments of inertia inertia[0] = tbonus[tri[i]].inertia[0]; inertia[1] = tbonus[tri[i]].inertia[1]; inertia[2] = tbonus[tri[i]].inertia[2]; // wbody = angular velocity in body frame MathExtra::quat_to_mat(tbonus[tri[i]].quat,rot); MathExtra::transpose_matvec(rot,angmom[i],wbody); wbody[0] /= inertia[0]; wbody[1] /= inertia[1]; wbody[2] /= inertia[2]; erotate += inertia[0]*wbody[0]*wbody[0] + inertia[1]*wbody[1]*wbody[1] + inertia[2]*wbody[2]*wbody[2]; } } MPI_Allreduce(&erotate,&scalar,1,MPI_DOUBLE,MPI_SUM,world); scalar *= pfactor; return scalar; }