2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
2007-01-30 08:22:05 +08:00
|
|
|
http://lammps.sandia.gov, Sandia National Laboratories
|
|
|
|
Steve Plimpton, sjplimp@sandia.gov
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
Copyright (2003) Sandia Corporation. Under the terms of Contract
|
|
|
|
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
|
2012-06-07 06:47:51 +08:00
|
|
|
certain rights in this software. This software is distributed under
|
2006-09-28 03:51:33 +08:00
|
|
|
the GNU General Public License.
|
|
|
|
|
|
|
|
See the README file in the top-level LAMMPS directory.
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
2015-10-31 04:04:06 +08:00
|
|
|
#include <mpi.h>
|
|
|
|
#include <math.h>
|
|
|
|
#include <stdio.h>
|
|
|
|
#include <string.h>
|
|
|
|
#include <stdlib.h>
|
2006-09-28 03:51:33 +08:00
|
|
|
#include "group.h"
|
|
|
|
#include "domain.h"
|
2007-01-30 08:22:05 +08:00
|
|
|
#include "atom.h"
|
2007-04-21 03:03:20 +08:00
|
|
|
#include "force.h"
|
2014-06-28 07:28:02 +08:00
|
|
|
#include "comm.h"
|
2006-09-28 03:51:33 +08:00
|
|
|
#include "region.h"
|
2009-03-17 23:40:57 +08:00
|
|
|
#include "modify.h"
|
|
|
|
#include "fix.h"
|
|
|
|
#include "compute.h"
|
|
|
|
#include "output.h"
|
2012-08-28 23:07:23 +08:00
|
|
|
#include "input.h"
|
|
|
|
#include "variable.h"
|
2009-03-17 23:40:57 +08:00
|
|
|
#include "dump.h"
|
2016-09-07 05:55:39 +08:00
|
|
|
#include "math_extra.h"
|
2012-08-28 23:07:23 +08:00
|
|
|
#include "memory.h"
|
2006-09-28 03:51:33 +08:00
|
|
|
#include "error.h"
|
2007-01-30 08:22:05 +08:00
|
|
|
|
2014-06-28 07:28:02 +08:00
|
|
|
#include <map>
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
using namespace LAMMPS_NS;
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
#define MAX_GROUP 32
|
2016-09-07 05:55:39 +08:00
|
|
|
#define EPSILON 1.0e-6
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2007-06-20 21:17:59 +08:00
|
|
|
enum{TYPE,MOLECULE,ID};
|
|
|
|
enum{LT,LE,GT,GE,EQ,NEQ,BETWEEN};
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
#define BIG 1.0e20
|
|
|
|
|
2014-06-28 07:28:02 +08:00
|
|
|
// allocate space for static class variable
|
|
|
|
|
|
|
|
Group *Group::cptr;
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
initialize group memory
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
Group::Group(LAMMPS *lmp) : Pointers(lmp)
|
2006-09-28 03:51:33 +08:00
|
|
|
{
|
|
|
|
MPI_Comm_rank(world,&me);
|
|
|
|
|
2009-03-17 23:40:57 +08:00
|
|
|
names = new char*[MAX_GROUP];
|
2006-09-28 03:51:33 +08:00
|
|
|
bitmask = new int[MAX_GROUP];
|
|
|
|
inversemask = new int[MAX_GROUP];
|
2014-02-06 01:30:23 +08:00
|
|
|
dynamic = new int[MAX_GROUP];
|
2009-03-17 23:40:57 +08:00
|
|
|
|
|
|
|
for (int i = 0; i < MAX_GROUP; i++) names[i] = NULL;
|
2006-09-28 03:51:33 +08:00
|
|
|
for (int i = 0; i < MAX_GROUP; i++) bitmask[i] = 1 << i;
|
|
|
|
for (int i = 0; i < MAX_GROUP; i++) inversemask[i] = bitmask[i] ^ ~0;
|
2014-02-06 01:30:23 +08:00
|
|
|
for (int i = 0; i < MAX_GROUP; i++) dynamic[i] = 0;
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
// create "all" group
|
|
|
|
|
2007-10-04 00:22:30 +08:00
|
|
|
char *str = (char *) "all";
|
2006-09-28 03:51:33 +08:00
|
|
|
int n = strlen(str) + 1;
|
2009-03-17 23:40:57 +08:00
|
|
|
names[0] = new char[n];
|
2006-09-28 03:51:33 +08:00
|
|
|
strcpy(names[0],str);
|
|
|
|
ngroup = 1;
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
free all memory
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
Group::~Group()
|
|
|
|
{
|
2009-03-17 23:40:57 +08:00
|
|
|
for (int i = 0; i < MAX_GROUP; i++) delete [] names[i];
|
|
|
|
delete [] names;
|
2006-09-28 03:51:33 +08:00
|
|
|
delete [] bitmask;
|
|
|
|
delete [] inversemask;
|
2014-02-06 01:30:23 +08:00
|
|
|
delete [] dynamic;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
assign atoms to a new or existing group
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::assign(int narg, char **arg)
|
|
|
|
{
|
|
|
|
int i;
|
|
|
|
|
2012-06-07 06:47:51 +08:00
|
|
|
if (domain->box_exist == 0)
|
2011-09-24 02:06:55 +08:00
|
|
|
error->all(FLERR,"Group command before simulation box is defined");
|
|
|
|
if (narg < 2) error->all(FLERR,"Illegal group command");
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2009-03-17 23:40:57 +08:00
|
|
|
// delete the group if not being used elsewhere
|
|
|
|
// clear mask of each atom assigned to this group
|
|
|
|
|
|
|
|
if (strcmp(arg[1],"delete") == 0) {
|
|
|
|
int igroup = find(arg[0]);
|
2011-09-24 02:06:55 +08:00
|
|
|
if (igroup == -1) error->all(FLERR,"Could not find group delete group ID");
|
|
|
|
if (igroup == 0) error->all(FLERR,"Cannot delete group all");
|
2009-03-17 23:40:57 +08:00
|
|
|
for (i = 0; i < modify->nfix; i++)
|
|
|
|
if (modify->fix[i]->igroup == igroup)
|
2012-06-07 06:47:51 +08:00
|
|
|
error->all(FLERR,"Cannot delete group currently used by a fix");
|
2009-03-17 23:40:57 +08:00
|
|
|
for (i = 0; i < modify->ncompute; i++)
|
|
|
|
if (modify->compute[i]->igroup == igroup)
|
2012-06-07 06:47:51 +08:00
|
|
|
error->all(FLERR,"Cannot delete group currently used by a compute");
|
2009-03-17 23:40:57 +08:00
|
|
|
for (i = 0; i < output->ndump; i++)
|
|
|
|
if (output->dump[i]->igroup == igroup)
|
2012-06-07 06:47:51 +08:00
|
|
|
error->all(FLERR,"Cannot delete group currently used by a dump");
|
2009-03-17 23:40:57 +08:00
|
|
|
if (atom->firstgroupname && strcmp(arg[0],atom->firstgroupname) == 0)
|
2012-02-16 08:35:51 +08:00
|
|
|
error->all(FLERR,
|
2012-06-07 06:47:51 +08:00
|
|
|
"Cannot delete group currently used by atom_modify first");
|
2009-03-17 23:40:57 +08:00
|
|
|
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
int bits = inversemask[igroup];
|
|
|
|
for (i = 0; i < nlocal; i++) mask[i] &= bits;
|
|
|
|
|
2014-02-06 01:30:23 +08:00
|
|
|
if (dynamic[igroup]) {
|
|
|
|
int n = strlen("GROUP_") + strlen(names[igroup]) + 1;
|
|
|
|
char *fixID = new char[n];
|
|
|
|
sprintf(fixID,"GROUP_%s",names[igroup]);
|
|
|
|
modify->delete_fix(fixID);
|
|
|
|
delete [] fixID;
|
|
|
|
}
|
|
|
|
|
2009-03-17 23:40:57 +08:00
|
|
|
delete [] names[igroup];
|
|
|
|
names[igroup] = NULL;
|
2014-02-06 01:30:23 +08:00
|
|
|
dynamic[igroup] = 0;
|
2009-03-17 23:40:57 +08:00
|
|
|
ngroup--;
|
|
|
|
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
2014-06-28 07:28:02 +08:00
|
|
|
// clear the group
|
|
|
|
|
|
|
|
if (strcmp(arg[1],"clear") == 0) {
|
|
|
|
int igroup = find(arg[0]);
|
|
|
|
if (igroup == -1) error->all (FLERR,"Could not find group clear group ID");
|
|
|
|
if (igroup == 0) error->all (FLERR,"Cannot clear group all");
|
|
|
|
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
int bits = inversemask[igroup];
|
|
|
|
for (i = 0; i < nlocal; i++) mask[i] &= bits;
|
|
|
|
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
// find group in existing list
|
2009-03-17 23:40:57 +08:00
|
|
|
// add a new group if igroup = -1
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
int igroup = find(arg[0]);
|
|
|
|
|
|
|
|
if (igroup == -1) {
|
2011-09-24 02:06:55 +08:00
|
|
|
if (ngroup == MAX_GROUP) error->all(FLERR,"Too many groups");
|
2009-03-17 23:40:57 +08:00
|
|
|
igroup = find_unused();
|
2006-09-28 03:51:33 +08:00
|
|
|
int n = strlen(arg[0]) + 1;
|
2009-03-17 23:40:57 +08:00
|
|
|
names[igroup] = new char[n];
|
2006-09-28 03:51:33 +08:00
|
|
|
strcpy(names[igroup],arg[0]);
|
2009-03-17 23:40:57 +08:00
|
|
|
ngroup++;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
int bit = bitmask[igroup];
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
// style = region
|
2007-06-20 21:17:59 +08:00
|
|
|
// add to group if atom is in region
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
if (strcmp(arg[1],"region") == 0) {
|
|
|
|
|
2011-09-24 02:06:55 +08:00
|
|
|
if (narg != 3) error->all(FLERR,"Illegal group command");
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2007-06-20 21:17:59 +08:00
|
|
|
int iregion = domain->find_region(arg[2]);
|
2011-09-24 02:06:55 +08:00
|
|
|
if (iregion == -1) error->all(FLERR,"Group region ID does not exist");
|
2013-08-02 22:50:05 +08:00
|
|
|
domain->regions[iregion]->init();
|
2014-05-02 22:28:56 +08:00
|
|
|
domain->regions[iregion]->prematch();
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (domain->regions[iregion]->match(x[i][0],x[i][1],x[i][2]))
|
2012-06-07 06:47:51 +08:00
|
|
|
mask[i] |= bit;
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2013-05-11 02:23:33 +08:00
|
|
|
// style = type, molecule, id
|
2014-02-06 01:30:23 +08:00
|
|
|
// add to group if atom matches type/molecule/id or condition
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"type") == 0 || strcmp(arg[1],"molecule") == 0 ||
|
2012-06-07 06:47:51 +08:00
|
|
|
strcmp(arg[1],"id") == 0) {
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2011-09-24 02:06:55 +08:00
|
|
|
if (narg < 3) error->all(FLERR,"Illegal group command");
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
int category;
|
|
|
|
if (strcmp(arg[1],"type") == 0) category = TYPE;
|
|
|
|
else if (strcmp(arg[1],"molecule") == 0) category = MOLECULE;
|
|
|
|
else if (strcmp(arg[1],"id") == 0) category = ID;
|
|
|
|
|
2013-05-11 02:23:33 +08:00
|
|
|
// args = logical condition
|
2015-10-31 04:04:06 +08:00
|
|
|
|
2013-05-11 02:23:33 +08:00
|
|
|
if (narg > 3 &&
|
|
|
|
(strcmp(arg[2],"<") == 0 || strcmp(arg[2],">") == 0 ||
|
|
|
|
strcmp(arg[2],"<=") == 0 || strcmp(arg[2],">=") == 0 ||
|
2014-06-02 22:24:56 +08:00
|
|
|
strcmp(arg[2],"==") == 0 || strcmp(arg[2],"!=") == 0 ||
|
2013-05-11 02:23:33 +08:00
|
|
|
strcmp(arg[2],"<>") == 0)) {
|
|
|
|
|
2014-04-11 22:46:49 +08:00
|
|
|
int condition = -1;
|
2013-05-11 02:23:33 +08:00
|
|
|
if (strcmp(arg[2],"<") == 0) condition = LT;
|
|
|
|
else if (strcmp(arg[2],"<=") == 0) condition = LE;
|
|
|
|
else if (strcmp(arg[2],">") == 0) condition = GT;
|
|
|
|
else if (strcmp(arg[2],">=") == 0) condition = GE;
|
|
|
|
else if (strcmp(arg[2],"==") == 0) condition = EQ;
|
|
|
|
else if (strcmp(arg[2],"!=") == 0) condition = NEQ;
|
|
|
|
else if (strcmp(arg[2],"<>") == 0) condition = BETWEEN;
|
|
|
|
else error->all(FLERR,"Illegal group command");
|
2015-10-31 04:04:06 +08:00
|
|
|
|
2014-01-18 02:43:09 +08:00
|
|
|
tagint bound1,bound2;
|
2014-04-30 06:43:28 +08:00
|
|
|
bound1 = force->tnumeric(FLERR,arg[3]);
|
2013-05-11 02:23:33 +08:00
|
|
|
bound2 = -1;
|
|
|
|
|
|
|
|
if (condition == BETWEEN) {
|
|
|
|
if (narg != 5) error->all(FLERR,"Illegal group command");
|
2014-04-30 06:43:28 +08:00
|
|
|
bound2 = force->tnumeric(FLERR,arg[4]);
|
2013-05-11 02:23:33 +08:00
|
|
|
} else if (narg != 4) error->all(FLERR,"Illegal group command");
|
|
|
|
|
2014-01-23 01:44:32 +08:00
|
|
|
int *attribute = NULL;
|
|
|
|
tagint *tattribute = NULL;
|
2013-05-11 02:23:33 +08:00
|
|
|
if (category == TYPE) attribute = atom->type;
|
2014-01-23 01:44:32 +08:00
|
|
|
else if (category == MOLECULE) tattribute = atom->molecule;
|
|
|
|
else if (category == ID) tattribute = atom->tag;
|
2013-05-11 02:23:33 +08:00
|
|
|
|
|
|
|
// add to group if meets condition
|
|
|
|
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute) {
|
|
|
|
if (condition == LT) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute[i] < bound1) mask[i] |= bit;
|
|
|
|
} else if (condition == LE) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute[i] <= bound1) mask[i] |= bit;
|
|
|
|
} else if (condition == GT) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute[i] > bound1) mask[i] |= bit;
|
|
|
|
} else if (condition == GE) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute[i] >= bound1) mask[i] |= bit;
|
|
|
|
} else if (condition == EQ) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute[i] == bound1) mask[i] |= bit;
|
|
|
|
} else if (condition == NEQ) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute[i] != bound1) mask[i] |= bit;
|
|
|
|
} else if (condition == BETWEEN) {
|
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (attribute[i] >= bound1 && attribute[i] <= bound2)
|
|
|
|
mask[i] |= bit;
|
|
|
|
}
|
|
|
|
} else {
|
|
|
|
if (condition == LT) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] < bound1) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
} else if (condition == LE) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] <= bound1) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
} else if (condition == GT) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] > bound1) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
} else if (condition == GE) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] >= bound1) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
} else if (condition == EQ) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] == bound1) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
} else if (condition == NEQ) {
|
2015-10-31 04:04:06 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] != bound1) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
} else if (condition == BETWEEN) {
|
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] >= bound1 && tattribute[i] <= bound2)
|
2014-01-18 02:43:09 +08:00
|
|
|
mask[i] |= bit;
|
|
|
|
}
|
2013-05-11 02:23:33 +08:00
|
|
|
}
|
2014-01-18 02:43:09 +08:00
|
|
|
|
2013-05-11 02:23:33 +08:00
|
|
|
// args = list of values
|
2015-10-31 04:04:06 +08:00
|
|
|
|
2013-05-11 02:23:33 +08:00
|
|
|
} else {
|
2014-01-23 01:44:32 +08:00
|
|
|
int *attribute = NULL;
|
|
|
|
tagint *tattribute = NULL;
|
2013-05-11 02:23:33 +08:00
|
|
|
if (category == TYPE) attribute = atom->type;
|
2014-01-23 01:44:32 +08:00
|
|
|
else if (category == MOLECULE) tattribute = atom->molecule;
|
|
|
|
else if (category == ID) tattribute = atom->tag;
|
2015-10-31 04:04:06 +08:00
|
|
|
|
2013-05-11 02:23:33 +08:00
|
|
|
char *ptr;
|
2014-01-18 02:43:09 +08:00
|
|
|
tagint start,stop,delta;
|
2013-05-11 02:23:33 +08:00
|
|
|
|
|
|
|
for (int iarg = 2; iarg < narg; iarg++) {
|
2015-09-03 04:39:55 +08:00
|
|
|
delta = 1;
|
2014-01-23 01:44:32 +08:00
|
|
|
if (strchr(arg[iarg],':')) {
|
2015-10-31 04:04:06 +08:00
|
|
|
ptr = strtok(arg[iarg],":");
|
|
|
|
start = force->tnumeric(FLERR,ptr);
|
|
|
|
ptr = strtok(NULL,":");
|
|
|
|
stop = force->tnumeric(FLERR,ptr);
|
2014-01-23 01:44:32 +08:00
|
|
|
ptr = strtok(NULL,":");
|
2014-04-30 06:43:28 +08:00
|
|
|
if (ptr) delta = force->tnumeric(FLERR,ptr);
|
2013-05-11 02:23:33 +08:00
|
|
|
} else {
|
2014-04-30 06:43:28 +08:00
|
|
|
start = stop = force->tnumeric(FLERR,arg[iarg]);
|
2013-05-11 02:23:33 +08:00
|
|
|
}
|
2015-09-03 04:39:55 +08:00
|
|
|
if (delta < 1)
|
|
|
|
error->all(FLERR,"Illegal range increment value");
|
2013-05-11 02:23:33 +08:00
|
|
|
|
|
|
|
// add to group if attribute matches value or sequence
|
2015-10-31 04:04:06 +08:00
|
|
|
|
2014-01-18 02:43:09 +08:00
|
|
|
if (attribute) {
|
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (attribute[i] >= start && attribute[i] <= stop &&
|
|
|
|
(attribute[i]-start) % delta == 0) mask[i] |= bit;
|
|
|
|
} else {
|
|
|
|
for (i = 0; i < nlocal; i++)
|
2014-01-23 01:44:32 +08:00
|
|
|
if (tattribute[i] >= start && tattribute[i] <= stop &&
|
|
|
|
(tattribute[i]-start) % delta == 0) mask[i] |= bit;
|
2014-01-18 02:43:09 +08:00
|
|
|
}
|
2013-05-11 02:23:33 +08:00
|
|
|
}
|
|
|
|
}
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2012-08-28 23:07:23 +08:00
|
|
|
// style = variable
|
2014-02-06 01:30:23 +08:00
|
|
|
// add to group if atom-atyle variable is non-zero
|
2012-08-28 23:07:23 +08:00
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"variable") == 0) {
|
|
|
|
|
|
|
|
int ivar = input->variable->find(arg[2]);
|
|
|
|
if (ivar < 0) error->all(FLERR,"Variable name for group does not exist");
|
|
|
|
if (!input->variable->atomstyle(ivar))
|
|
|
|
error->all(FLERR,"Variable for group is invalid style");
|
|
|
|
|
|
|
|
double *aflag;
|
2015-10-31 04:04:06 +08:00
|
|
|
|
2012-08-28 23:07:23 +08:00
|
|
|
// aflag = evaluation of per-atom variable
|
|
|
|
|
|
|
|
memory->create(aflag,nlocal,"group:aflag");
|
|
|
|
input->variable->compute_atom(ivar,0,aflag,1,0);
|
|
|
|
|
|
|
|
// add to group if per-atom variable evaluated to non-zero
|
|
|
|
|
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (aflag[i] != 0.0) mask[i] |= bit;
|
|
|
|
|
|
|
|
memory->destroy(aflag);
|
|
|
|
|
2014-06-28 07:28:02 +08:00
|
|
|
// style = include
|
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"include") == 0) {
|
|
|
|
|
|
|
|
if (narg != 3) error->all(FLERR,"Illegal group command");
|
2015-10-31 04:04:06 +08:00
|
|
|
if (strcmp(arg[2],"molecule") != 0)
|
2014-06-28 07:28:02 +08:00
|
|
|
error->all(FLERR,"Illegal group command");
|
|
|
|
|
|
|
|
add_molecules(igroup,bit);
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
// style = subtract
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"subtract") == 0) {
|
|
|
|
|
2011-09-24 02:06:55 +08:00
|
|
|
if (narg < 4) error->all(FLERR,"Illegal group command");
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
int length = narg-2;
|
|
|
|
int *list = new int[length];
|
|
|
|
|
|
|
|
int jgroup;
|
|
|
|
for (int iarg = 2; iarg < narg; iarg++) {
|
|
|
|
jgroup = find(arg[iarg]);
|
2011-09-24 02:06:55 +08:00
|
|
|
if (jgroup == -1) error->all(FLERR,"Group ID does not exist");
|
2015-10-31 04:04:06 +08:00
|
|
|
if (dynamic[jgroup])
|
2014-02-06 01:30:23 +08:00
|
|
|
error->all(FLERR,"Cannot subtract groups using a dynamic group");
|
2006-09-28 03:51:33 +08:00
|
|
|
list[iarg-2] = jgroup;
|
|
|
|
}
|
|
|
|
|
|
|
|
// add to group if in 1st group in list
|
|
|
|
|
|
|
|
int otherbit = bitmask[list[0]];
|
|
|
|
|
2012-06-07 06:47:51 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2006-09-28 03:51:33 +08:00
|
|
|
if (mask[i] & otherbit) mask[i] |= bit;
|
|
|
|
|
|
|
|
// remove atoms if they are in any of the other groups
|
|
|
|
// AND with inverse mask removes the atom from group
|
|
|
|
|
|
|
|
int inverse = inversemask[igroup];
|
|
|
|
|
|
|
|
for (int ilist = 1; ilist < length; ilist++) {
|
|
|
|
otherbit = bitmask[list[ilist]];
|
2012-06-07 06:47:51 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & otherbit) mask[i] &= inverse;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
delete [] list;
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
// style = union
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"union") == 0) {
|
|
|
|
|
2011-09-24 02:06:55 +08:00
|
|
|
if (narg < 3) error->all(FLERR,"Illegal group command");
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
int length = narg-2;
|
|
|
|
int *list = new int[length];
|
|
|
|
|
|
|
|
int jgroup;
|
|
|
|
for (int iarg = 2; iarg < narg; iarg++) {
|
|
|
|
jgroup = find(arg[iarg]);
|
2011-09-24 02:06:55 +08:00
|
|
|
if (jgroup == -1) error->all(FLERR,"Group ID does not exist");
|
2015-10-31 04:04:06 +08:00
|
|
|
if (dynamic[jgroup])
|
2014-02-06 01:30:23 +08:00
|
|
|
error->all(FLERR,"Cannot union groups using a dynamic group");
|
2006-09-28 03:51:33 +08:00
|
|
|
list[iarg-2] = jgroup;
|
|
|
|
}
|
|
|
|
|
|
|
|
// add to group if in any other group in list
|
|
|
|
|
|
|
|
int otherbit;
|
|
|
|
|
|
|
|
for (int ilist = 0; ilist < length; ilist++) {
|
|
|
|
otherbit = bitmask[list[ilist]];
|
2012-06-07 06:47:51 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & otherbit) mask[i] |= bit;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
delete [] list;
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
// style = intersect
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"intersect") == 0) {
|
|
|
|
|
2011-09-24 02:06:55 +08:00
|
|
|
if (narg < 4) error->all(FLERR,"Illegal group command");
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
int length = narg-2;
|
|
|
|
int *list = new int[length];
|
|
|
|
|
|
|
|
int jgroup;
|
|
|
|
for (int iarg = 2; iarg < narg; iarg++) {
|
|
|
|
jgroup = find(arg[iarg]);
|
2011-09-24 02:06:55 +08:00
|
|
|
if (jgroup == -1) error->all(FLERR,"Group ID does not exist");
|
2015-10-31 04:04:06 +08:00
|
|
|
if (dynamic[jgroup])
|
2014-02-06 01:30:23 +08:00
|
|
|
error->all(FLERR,"Cannot intersect groups using a dynamic group");
|
2006-09-28 03:51:33 +08:00
|
|
|
list[iarg-2] = jgroup;
|
|
|
|
}
|
|
|
|
|
|
|
|
// add to group if in all groups in list
|
|
|
|
|
|
|
|
int otherbit,ok,ilist;
|
|
|
|
|
|
|
|
for (i = 0; i < nlocal; i++) {
|
|
|
|
ok = 1;
|
|
|
|
for (ilist = 0; ilist < length; ilist++) {
|
2012-06-07 06:47:51 +08:00
|
|
|
otherbit = bitmask[list[ilist]];
|
|
|
|
if ((mask[i] & otherbit) == 0) ok = 0;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
if (ok) mask[i] |= bit;
|
|
|
|
}
|
|
|
|
|
|
|
|
delete [] list;
|
|
|
|
|
2014-02-06 01:30:23 +08:00
|
|
|
// style = dynamic
|
|
|
|
// create a new FixGroup to dynamically determine atoms in group
|
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"dynamic") == 0) {
|
|
|
|
|
|
|
|
if (narg < 4) error->all(FLERR,"Illegal group command");
|
2015-10-31 04:04:06 +08:00
|
|
|
if (strcmp(arg[0],arg[2]) == 0)
|
2014-02-06 03:38:49 +08:00
|
|
|
error->all(FLERR,"Group dynamic cannot reference itself");
|
2015-10-31 04:04:06 +08:00
|
|
|
if (find(arg[2]) < 0)
|
2014-02-06 03:38:49 +08:00
|
|
|
error->all(FLERR,"Group dynamic parent group does not exist");
|
2014-02-06 05:31:26 +08:00
|
|
|
if (igroup == 0) error->all(FLERR,"Group all cannot be made dynamic");
|
2014-02-06 01:30:23 +08:00
|
|
|
|
2014-02-06 03:38:49 +08:00
|
|
|
// if group is already dynamic, delete existing FixGroup
|
2014-02-06 01:30:23 +08:00
|
|
|
|
|
|
|
if (dynamic[igroup]) {
|
|
|
|
int n = strlen("GROUP_") + strlen(names[igroup]) + 1;
|
|
|
|
char *fixID = new char[n];
|
|
|
|
sprintf(fixID,"GROUP_%s",names[igroup]);
|
|
|
|
modify->delete_fix(fixID);
|
|
|
|
delete [] fixID;
|
|
|
|
}
|
|
|
|
|
|
|
|
dynamic[igroup] = 1;
|
|
|
|
|
|
|
|
int n = strlen("GROUP_") + strlen(names[igroup]) + 1;
|
|
|
|
char *fixID = new char[n];
|
|
|
|
sprintf(fixID,"GROUP_%s",names[igroup]);
|
|
|
|
|
|
|
|
char **newarg = new char*[narg];
|
|
|
|
newarg[0] = fixID;
|
|
|
|
newarg[1] = arg[2];
|
|
|
|
newarg[2] = (char *) "GROUP";
|
|
|
|
for (int i = 3; i < narg; i++) newarg[i] = arg[i];
|
|
|
|
modify->add_fix(narg,newarg);
|
|
|
|
delete [] newarg;
|
|
|
|
delete [] fixID;
|
|
|
|
|
|
|
|
// style = static
|
|
|
|
// remove dynamic FixGroup if necessary
|
|
|
|
|
|
|
|
} else if (strcmp(arg[1],"static") == 0) {
|
|
|
|
|
|
|
|
if (narg != 2) error->all(FLERR,"Illegal group command");
|
|
|
|
|
|
|
|
if (dynamic[igroup]) {
|
|
|
|
int n = strlen("GROUP_") + strlen(names[igroup]) + 1;
|
|
|
|
char *fixID = new char[n];
|
|
|
|
sprintf(fixID,"GROUP_%s",names[igroup]);
|
|
|
|
modify->delete_fix(fixID);
|
|
|
|
delete [] fixID;
|
|
|
|
}
|
|
|
|
|
|
|
|
dynamic[igroup] = 0;
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
// not a valid group style
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2011-09-24 02:06:55 +08:00
|
|
|
} else error->all(FLERR,"Illegal group command");
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
// print stats for changed group
|
|
|
|
|
2008-06-07 01:31:38 +08:00
|
|
|
int n;
|
2006-09-28 03:51:33 +08:00
|
|
|
n = 0;
|
|
|
|
for (i = 0; i < nlocal; i++) if (mask[i] & bit) n++;
|
2008-06-07 01:31:38 +08:00
|
|
|
|
|
|
|
double rlocal = n;
|
|
|
|
double all;
|
|
|
|
MPI_Allreduce(&rlocal,&all,1,MPI_DOUBLE,MPI_SUM,world);
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
if (me == 0) {
|
2014-04-05 00:36:51 +08:00
|
|
|
if (dynamic[igroup]) {
|
2014-05-02 23:49:01 +08:00
|
|
|
if (screen) fprintf(screen,"dynamic group %s defined\n",names[igroup]);
|
|
|
|
if (logfile) fprintf(logfile,"dynamic group %s defined\n",names[igroup]);
|
2014-04-05 00:36:51 +08:00
|
|
|
} else {
|
2015-10-31 04:04:06 +08:00
|
|
|
if (screen)
|
2014-04-05 00:36:51 +08:00
|
|
|
fprintf(screen,"%.15g atoms in group %s\n",all,names[igroup]);
|
|
|
|
if (logfile)
|
|
|
|
fprintf(logfile,"%.15g atoms in group %s\n",all,names[igroup]);
|
|
|
|
}
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
}
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
add flagged atoms to a new or existing group
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::create(char *name, int *flag)
|
|
|
|
{
|
|
|
|
int i;
|
|
|
|
|
|
|
|
// find group in existing list
|
2009-03-17 23:40:57 +08:00
|
|
|
// add a new group if igroup = -1
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
int igroup = find(name);
|
|
|
|
|
|
|
|
if (igroup == -1) {
|
2011-09-24 02:06:55 +08:00
|
|
|
if (ngroup == MAX_GROUP) error->all(FLERR,"Too many groups");
|
2009-03-17 23:40:57 +08:00
|
|
|
igroup = find_unused();
|
2006-09-28 03:51:33 +08:00
|
|
|
int n = strlen(name) + 1;
|
2009-03-17 23:40:57 +08:00
|
|
|
names[igroup] = new char[n];
|
2006-09-28 03:51:33 +08:00
|
|
|
strcpy(names[igroup],name);
|
2009-03-17 23:40:57 +08:00
|
|
|
ngroup++;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
2009-03-17 23:40:57 +08:00
|
|
|
// add atoms to group whose flags are set
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
int bit = bitmask[igroup];
|
|
|
|
|
2012-06-07 06:47:51 +08:00
|
|
|
for (i = 0; i < nlocal; i++)
|
2006-09-28 03:51:33 +08:00
|
|
|
if (flag[i]) mask[i] |= bit;
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
return group index if name matches existing group, -1 if no such group
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
2007-10-04 00:22:30 +08:00
|
|
|
int Group::find(const char *name)
|
2006-09-28 03:51:33 +08:00
|
|
|
{
|
2009-03-17 23:40:57 +08:00
|
|
|
for (int igroup = 0; igroup < MAX_GROUP; igroup++)
|
|
|
|
if (names[igroup] && strcmp(name,names[igroup]) == 0) return igroup;
|
|
|
|
return -1;
|
|
|
|
}
|
|
|
|
|
2015-07-18 07:20:48 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
find group with name or create group if it doesn't exist
|
|
|
|
return group index
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
int Group::find_or_create(const char *name)
|
|
|
|
{
|
|
|
|
int igroup = find(name);
|
|
|
|
if (igroup >= 0) return igroup;
|
|
|
|
|
|
|
|
if (ngroup == MAX_GROUP) error->all(FLERR,"Too many groups");
|
|
|
|
igroup = find_unused();
|
|
|
|
int n = strlen(name) + 1;
|
|
|
|
names[igroup] = new char[n];
|
|
|
|
strcpy(names[igroup],name);
|
|
|
|
ngroup++;
|
|
|
|
|
|
|
|
return igroup;
|
|
|
|
}
|
|
|
|
|
2009-03-17 23:40:57 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
return index of first available group
|
|
|
|
should never be called when group limit has been reached
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
int Group::find_unused()
|
|
|
|
{
|
|
|
|
for (int igroup = 0; igroup < MAX_GROUP; igroup++)
|
|
|
|
if (names[igroup] == NULL) return igroup;
|
2006-09-28 03:51:33 +08:00
|
|
|
return -1;
|
|
|
|
}
|
|
|
|
|
2014-06-28 07:28:02 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
add atoms to group that are in same molecules as atoms already in group
|
|
|
|
do not include molID = 0
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::add_molecules(int igroup, int bit)
|
|
|
|
{
|
|
|
|
// hash = unique molecule IDs of atoms already in group
|
|
|
|
|
|
|
|
hash = new std::map<tagint,int>();
|
|
|
|
|
|
|
|
tagint *molecule = atom->molecule;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & bit) {
|
|
|
|
if (molecule[i] == 0) continue;
|
|
|
|
if (hash->find(molecule[i]) == hash->end()) (*hash)[molecule[i]] = 1;
|
|
|
|
}
|
|
|
|
|
|
|
|
// list = set of unique molecule IDs for atoms to add
|
|
|
|
// pass list to all other procs via comm->ring()
|
|
|
|
|
|
|
|
int n = hash->size();
|
|
|
|
tagint *list;
|
|
|
|
memory->create(list,n,"group:list");
|
|
|
|
|
|
|
|
n = 0;
|
|
|
|
std::map<tagint,int>::iterator pos;
|
|
|
|
for (pos = hash->begin(); pos != hash->end(); ++pos) list[n++] = pos->first;
|
|
|
|
|
|
|
|
cptr = this;
|
|
|
|
molbit = bit;
|
|
|
|
comm->ring(n,sizeof(tagint),list,1,molring,NULL);
|
|
|
|
|
|
|
|
delete hash;
|
|
|
|
memory->destroy(list);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
callback from comm->ring()
|
|
|
|
cbuf = list of N molecule IDs, put them in hash
|
|
|
|
loop over my atoms, if matches molecule ID in hash,
|
|
|
|
add atom to group flagged by molbit
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::molring(int n, char *cbuf)
|
|
|
|
{
|
|
|
|
tagint *list = (tagint *) cbuf;
|
|
|
|
std::map<tagint,int> *hash = cptr->hash;
|
|
|
|
int nlocal = cptr->atom->nlocal;
|
|
|
|
tagint *molecule = cptr->atom->molecule;
|
|
|
|
int *mask = cptr->atom->mask;
|
|
|
|
int molbit = cptr->molbit;
|
|
|
|
|
|
|
|
hash->clear();
|
|
|
|
for (int i = 0; i < n; i++) (*hash)[list[i]] = 1;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (hash->find(molecule[i]) != hash->end()) mask[i] |= molbit;
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
write group info to a restart file
|
|
|
|
only called by proc 0
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::write_restart(FILE *fp)
|
|
|
|
{
|
|
|
|
fwrite(&ngroup,sizeof(int),1,fp);
|
|
|
|
|
2009-03-18 22:36:01 +08:00
|
|
|
// use count to not change restart format with deleted groups
|
|
|
|
// remove this on next major release
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
int n;
|
2009-03-18 22:36:01 +08:00
|
|
|
int count = 0;
|
2009-03-17 23:40:57 +08:00
|
|
|
for (int i = 0; i < MAX_GROUP; i++) {
|
|
|
|
if (names[i]) n = strlen(names[i]) + 1;
|
|
|
|
else n = 0;
|
2006-09-28 03:51:33 +08:00
|
|
|
fwrite(&n,sizeof(int),1,fp);
|
2009-03-18 22:36:01 +08:00
|
|
|
if (n) {
|
|
|
|
fwrite(names[i],sizeof(char),n,fp);
|
|
|
|
count++;
|
|
|
|
}
|
|
|
|
if (count == ngroup) break;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
read group info from a restart file
|
|
|
|
proc 0 reads, bcast to all procs
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::read_restart(FILE *fp)
|
|
|
|
{
|
|
|
|
int i,n;
|
|
|
|
|
2009-03-17 23:40:57 +08:00
|
|
|
// delete existing group names
|
|
|
|
// atom masks will be overwritten by reading of restart file
|
|
|
|
|
|
|
|
for (i = 0; i < MAX_GROUP; i++) delete [] names[i];
|
2006-09-28 03:51:33 +08:00
|
|
|
|
|
|
|
if (me == 0) fread(&ngroup,sizeof(int),1,fp);
|
|
|
|
MPI_Bcast(&ngroup,1,MPI_INT,0,world);
|
2009-03-18 22:36:01 +08:00
|
|
|
|
|
|
|
// use count to not change restart format with deleted groups
|
|
|
|
// remove this on next major release
|
|
|
|
|
|
|
|
int count = 0;
|
2009-03-17 23:40:57 +08:00
|
|
|
for (i = 0; i < MAX_GROUP; i++) {
|
2009-03-18 22:36:01 +08:00
|
|
|
if (count == ngroup) {
|
|
|
|
names[i] = NULL;
|
|
|
|
continue;
|
|
|
|
}
|
2006-09-28 03:51:33 +08:00
|
|
|
if (me == 0) fread(&n,sizeof(int),1,fp);
|
|
|
|
MPI_Bcast(&n,1,MPI_INT,0,world);
|
2009-03-17 23:40:57 +08:00
|
|
|
if (n) {
|
|
|
|
names[i] = new char[n];
|
|
|
|
if (me == 0) fread(names[i],sizeof(char),n,fp);
|
|
|
|
MPI_Bcast(names[i],n,MPI_CHAR,0,world);
|
2009-03-18 22:36:01 +08:00
|
|
|
count++;
|
2009-03-17 23:40:57 +08:00
|
|
|
} else names[i] = NULL;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// ----------------------------------------------------------------------
|
|
|
|
// computations on a group of atoms
|
|
|
|
// ----------------------------------------------------------------------
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
count atoms in group
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
2011-01-05 05:08:27 +08:00
|
|
|
bigint Group::count(int igroup)
|
2006-09-28 03:51:33 +08:00
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
int n = 0;
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) n++;
|
|
|
|
|
2011-01-05 05:08:27 +08:00
|
|
|
bigint nsingle = n;
|
|
|
|
bigint nall;
|
2011-01-07 07:23:37 +08:00
|
|
|
MPI_Allreduce(&nsingle,&nall,1,MPI_LMP_BIGINT,MPI_SUM,world);
|
2006-09-28 03:51:33 +08:00
|
|
|
return nall;
|
|
|
|
}
|
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
count atoms in group and region
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
2011-01-05 05:08:27 +08:00
|
|
|
bigint Group::count(int igroup, int iregion)
|
2009-04-29 03:46:52 +08:00
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
int n = 0;
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) n++;
|
|
|
|
|
2011-01-05 05:08:27 +08:00
|
|
|
bigint nsingle = n;
|
|
|
|
bigint nall;
|
2011-01-07 07:23:37 +08:00
|
|
|
MPI_Allreduce(&nsingle,&nall,1,MPI_LMP_BIGINT,MPI_SUM,world);
|
2009-04-29 03:46:52 +08:00
|
|
|
return nall;
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the total mass of group of atoms
|
|
|
|
use either per-type mass or per-atom rmass
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::mass(int igroup)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
2009-04-29 03:46:52 +08:00
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2006-09-28 03:51:33 +08:00
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
2007-04-21 03:03:20 +08:00
|
|
|
double one = 0.0;
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) {
|
2006-09-28 03:51:33 +08:00
|
|
|
for (int i = 0; i < nlocal; i++)
|
2009-02-13 05:36:52 +08:00
|
|
|
if (mask[i] & groupbit) one += rmass[i];
|
2006-09-28 03:51:33 +08:00
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
2009-02-13 05:36:52 +08:00
|
|
|
if (mask[i] & groupbit) one += mass[type[i]];
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
2007-04-21 03:03:20 +08:00
|
|
|
double all;
|
|
|
|
MPI_Allreduce(&one,&all,1,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
return all;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the total mass of group of atoms in region
|
|
|
|
use either per-type mass or per-atom rmass
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::mass(int igroup, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double one = 0.0;
|
|
|
|
|
|
|
|
if (rmass) {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2]))
|
2012-06-07 06:47:51 +08:00
|
|
|
one += rmass[i];
|
2009-04-29 03:46:52 +08:00
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2]))
|
2012-06-07 06:47:51 +08:00
|
|
|
one += mass[type[i]];
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
double all;
|
|
|
|
MPI_Allreduce(&one,&all,1,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
return all;
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the total charge of group of atoms
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::charge(int igroup)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double *q = atom->q;
|
2009-04-29 03:46:52 +08:00
|
|
|
int *mask = atom->mask;
|
2006-09-28 03:51:33 +08:00
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double qone = 0.0;
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) qone += q[i];
|
|
|
|
|
|
|
|
double qall;
|
|
|
|
MPI_Allreduce(&qone,&qall,1,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
return qall;
|
|
|
|
}
|
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the total charge of group of atoms in region
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::charge(int igroup, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double *q = atom->q;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double qone = 0.0;
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2]))
|
|
|
|
qone += q[i];
|
|
|
|
|
|
|
|
double qall;
|
|
|
|
MPI_Allreduce(&qone,&qall,1,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
return qall;
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the coordinate bounds of the group of atoms
|
|
|
|
periodic images are not considered, so atoms are NOT unwrapped
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::bounds(int igroup, double *minmax)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double extent[6];
|
|
|
|
extent[0] = extent[2] = extent[4] = BIG;
|
|
|
|
extent[1] = extent[3] = extent[5] = -BIG;
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++) {
|
|
|
|
if (mask[i] & groupbit) {
|
|
|
|
extent[0] = MIN(extent[0],x[i][0]);
|
|
|
|
extent[1] = MAX(extent[1],x[i][0]);
|
|
|
|
extent[2] = MIN(extent[2],x[i][1]);
|
|
|
|
extent[3] = MAX(extent[3],x[i][1]);
|
|
|
|
extent[4] = MIN(extent[4],x[i][2]);
|
|
|
|
extent[5] = MAX(extent[5],x[i][2]);
|
|
|
|
}
|
|
|
|
}
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
// compute extent across all procs
|
|
|
|
// flip sign of MIN to do it in one Allreduce MAX
|
|
|
|
// set box by extent in shrink-wrapped dims
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
extent[0] = -extent[0];
|
|
|
|
extent[2] = -extent[2];
|
|
|
|
extent[4] = -extent[4];
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
MPI_Allreduce(extent,minmax,6,MPI_DOUBLE,MPI_MAX,world);
|
|
|
|
|
2007-06-20 21:17:59 +08:00
|
|
|
minmax[0] = -minmax[0];
|
|
|
|
minmax[2] = -minmax[2];
|
|
|
|
minmax[4] = -minmax[4];
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
2009-04-29 03:46:52 +08:00
|
|
|
compute the coordinate bounds of the group of atoms in region
|
|
|
|
periodic images are not considered, so atoms are NOT unwrapped
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::bounds(int igroup, double *minmax, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double extent[6];
|
|
|
|
extent[0] = extent[2] = extent[4] = BIG;
|
|
|
|
extent[1] = extent[3] = extent[5] = -BIG;
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++) {
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
|
|
|
extent[0] = MIN(extent[0],x[i][0]);
|
|
|
|
extent[1] = MAX(extent[1],x[i][0]);
|
|
|
|
extent[2] = MIN(extent[2],x[i][1]);
|
|
|
|
extent[3] = MAX(extent[3],x[i][1]);
|
|
|
|
extent[4] = MIN(extent[4],x[i][2]);
|
|
|
|
extent[5] = MAX(extent[5],x[i][2]);
|
|
|
|
}
|
|
|
|
}
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
// compute extent across all procs
|
|
|
|
// flip sign of MIN to do it in one Allreduce MAX
|
|
|
|
// set box by extent in shrink-wrapped dims
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
extent[0] = -extent[0];
|
|
|
|
extent[2] = -extent[2];
|
|
|
|
extent[4] = -extent[4];
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
MPI_Allreduce(extent,minmax,6,MPI_DOUBLE,MPI_MAX,world);
|
|
|
|
|
|
|
|
minmax[0] = -minmax[0];
|
|
|
|
minmax[2] = -minmax[2];
|
|
|
|
minmax[4] = -minmax[4];
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the center-of-mass coords of group of atoms
|
|
|
|
masstotal = total mass
|
2006-09-28 03:51:33 +08:00
|
|
|
return center-of-mass coords in cm[]
|
|
|
|
must unwrap atoms to compute center-of-mass correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::xcm(int igroup, double masstotal, double *cm)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2006-09-28 03:51:33 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double cmone[3];
|
|
|
|
cmone[0] = cmone[1] = cmone[2] = 0.0;
|
|
|
|
|
|
|
|
double massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
2006-09-28 03:51:33 +08:00
|
|
|
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) {
|
2006-09-28 03:51:33 +08:00
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = rmass[i];
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
cmone[0] += unwrap[0] * massone;
|
|
|
|
cmone[1] += unwrap[1] * massone;
|
|
|
|
cmone[2] += unwrap[2] * massone;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = mass[type[i]];
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
cmone[0] += unwrap[0] * massone;
|
|
|
|
cmone[1] += unwrap[1] * massone;
|
|
|
|
cmone[2] += unwrap[2] * massone;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(cmone,cm,3,MPI_DOUBLE,MPI_SUM,world);
|
2009-04-30 00:54:06 +08:00
|
|
|
if (masstotal > 0.0) {
|
|
|
|
cm[0] /= masstotal;
|
|
|
|
cm[1] /= masstotal;
|
|
|
|
cm[2] /= masstotal;
|
|
|
|
}
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
2009-04-29 03:46:52 +08:00
|
|
|
compute the center-of-mass coords of group of atoms in region
|
|
|
|
mastotal = total mass
|
|
|
|
return center-of-mass coords in cm[]
|
|
|
|
must unwrap atoms to compute center-of-mass correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::xcm(int igroup, double masstotal, double *cm, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2009-04-29 03:46:52 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double cmone[3];
|
|
|
|
cmone[0] = cmone[1] = cmone[2] = 0.0;
|
|
|
|
|
|
|
|
double massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
if (rmass) {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = rmass[i];
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
cmone[0] += unwrap[0] * massone;
|
|
|
|
cmone[1] += unwrap[1] * massone;
|
|
|
|
cmone[2] += unwrap[2] * massone;
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = mass[type[i]];
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
cmone[0] += unwrap[0] * massone;
|
|
|
|
cmone[1] += unwrap[1] * massone;
|
|
|
|
cmone[2] += unwrap[2] * massone;
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(cmone,cm,3,MPI_DOUBLE,MPI_SUM,world);
|
2009-04-30 00:54:06 +08:00
|
|
|
if (masstotal > 0.0) {
|
|
|
|
cm[0] /= masstotal;
|
|
|
|
cm[1] /= masstotal;
|
|
|
|
cm[2] /= masstotal;
|
|
|
|
}
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the center-of-mass velocity of group of atoms
|
|
|
|
masstotal = total mass
|
2006-09-28 03:51:33 +08:00
|
|
|
return center-of-mass velocity in cm[]
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::vcm(int igroup, double masstotal, double *cm)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **v = atom->v;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double p[3],massone;
|
|
|
|
p[0] = p[1] = p[2] = 0.0;
|
|
|
|
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) {
|
2006-09-28 03:51:33 +08:00
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = rmass[i];
|
|
|
|
p[0] += v[i][0]*massone;
|
|
|
|
p[1] += v[i][1]*massone;
|
|
|
|
p[2] += v[i][2]*massone;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = mass[type[i]];
|
|
|
|
p[0] += v[i][0]*massone;
|
|
|
|
p[1] += v[i][1]*massone;
|
|
|
|
p[2] += v[i][2]*massone;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(p,cm,3,MPI_DOUBLE,MPI_SUM,world);
|
2009-04-30 00:54:06 +08:00
|
|
|
if (masstotal > 0.0) {
|
|
|
|
cm[0] /= masstotal;
|
|
|
|
cm[1] /= masstotal;
|
|
|
|
cm[2] /= masstotal;
|
|
|
|
}
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the center-of-mass velocity of group of atoms in region
|
|
|
|
masstotal = total mass
|
|
|
|
return center-of-mass velocity in cm[]
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::vcm(int igroup, double masstotal, double *cm, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **v = atom->v;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double p[3],massone;
|
|
|
|
p[0] = p[1] = p[2] = 0.0;
|
|
|
|
|
|
|
|
if (rmass) {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = rmass[i];
|
|
|
|
p[0] += v[i][0]*massone;
|
|
|
|
p[1] += v[i][1]*massone;
|
|
|
|
p[2] += v[i][2]*massone;
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-06-07 06:47:51 +08:00
|
|
|
massone = mass[type[i]];
|
|
|
|
p[0] += v[i][0]*massone;
|
|
|
|
p[1] += v[i][1]*massone;
|
|
|
|
p[2] += v[i][2]*massone;
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(p,cm,3,MPI_DOUBLE,MPI_SUM,world);
|
2009-04-30 00:54:06 +08:00
|
|
|
if (masstotal > 0.0) {
|
|
|
|
cm[0] /= masstotal;
|
|
|
|
cm[1] /= masstotal;
|
|
|
|
cm[2] /= masstotal;
|
|
|
|
}
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
|
2007-03-23 01:05:51 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the total force on group of atoms
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::fcm(int igroup, double *cm)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **f = atom->f;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double flocal[3];
|
|
|
|
flocal[0] = flocal[1] = flocal[2] = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
|
|
|
flocal[0] += f[i][0];
|
|
|
|
flocal[1] += f[i][1];
|
|
|
|
flocal[2] += f[i][2];
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(flocal,cm,3,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
2007-04-21 03:03:20 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
2009-04-29 03:46:52 +08:00
|
|
|
compute the total force on group of atoms in region
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::fcm(int igroup, double *cm, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **f = atom->f;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double flocal[3];
|
|
|
|
flocal[0] = flocal[1] = flocal[2] = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
|
|
|
flocal[0] += f[i][0];
|
|
|
|
flocal[1] += f[i][1];
|
|
|
|
flocal[2] += f[i][2];
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(flocal,cm,3,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the total kinetic energy of group of atoms and return it
|
2007-04-21 03:03:20 +08:00
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::ke(int igroup)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **v = atom->v;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double one = 0.0;
|
|
|
|
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) {
|
2007-04-21 03:03:20 +08:00
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit)
|
2012-06-07 06:47:51 +08:00
|
|
|
one += (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]) *
|
|
|
|
rmass[i];
|
2007-04-21 03:03:20 +08:00
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit)
|
2012-06-07 06:47:51 +08:00
|
|
|
one += (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]) *
|
|
|
|
mass[type[i]];
|
2007-04-21 03:03:20 +08:00
|
|
|
}
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2007-04-21 03:03:20 +08:00
|
|
|
double all;
|
|
|
|
MPI_Allreduce(&one,&all,1,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
all *= 0.5 * force->mvv2e;
|
|
|
|
return all;
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
2009-04-29 03:46:52 +08:00
|
|
|
compute the total kinetic energy of group of atoms in region and return it
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::ke(int igroup, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **v = atom->v;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double one = 0.0;
|
|
|
|
|
|
|
|
if (rmass) {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2]))
|
2012-06-07 06:47:51 +08:00
|
|
|
one += (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]) *
|
|
|
|
rmass[i];
|
2009-04-29 03:46:52 +08:00
|
|
|
} else {
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2]))
|
2012-06-07 06:47:51 +08:00
|
|
|
one += (v[i][0]*v[i][0] + v[i][1]*v[i][1] + v[i][2]*v[i][2]) *
|
|
|
|
mass[type[i]];
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-29 03:46:52 +08:00
|
|
|
double all;
|
|
|
|
MPI_Allreduce(&one,&all,1,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
all *= 0.5 * force->mvv2e;
|
|
|
|
return all;
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the radius-of-gyration of group of atoms
|
|
|
|
around center-of-mass cm
|
2006-09-28 03:51:33 +08:00
|
|
|
must unwrap atoms to compute Rg correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::gyration(int igroup, double masstotal, double *cm)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2006-09-28 03:51:33 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz,massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
2006-09-28 03:51:33 +08:00
|
|
|
double rg = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2009-04-29 03:46:52 +08:00
|
|
|
if (rmass) massone = rmass[i];
|
|
|
|
else massone = mass[type[i]];
|
|
|
|
rg += (dx*dx + dy*dy + dz*dz) * massone;
|
|
|
|
}
|
|
|
|
double rg_all;
|
|
|
|
MPI_Allreduce(&rg,&rg_all,1,MPI_DOUBLE,MPI_SUM,world);
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-30 00:54:06 +08:00
|
|
|
if (masstotal > 0.0) return sqrt(rg_all/masstotal);
|
|
|
|
return 0.0;
|
2009-04-29 03:46:52 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the radius-of-gyration of group of atoms in region
|
|
|
|
around center-of-mass cm
|
|
|
|
must unwrap atoms to compute Rg correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
double Group::gyration(int igroup, double masstotal, double *cm, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2009-04-29 03:46:52 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2009-04-29 03:46:52 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz,massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
2009-04-29 03:46:52 +08:00
|
|
|
double rg = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) massone = rmass[i];
|
|
|
|
else massone = mass[type[i]];
|
2006-09-28 03:51:33 +08:00
|
|
|
rg += (dx*dx + dy*dy + dz*dz) * massone;
|
|
|
|
}
|
|
|
|
double rg_all;
|
|
|
|
MPI_Allreduce(&rg,&rg_all,1,MPI_DOUBLE,MPI_SUM,world);
|
2012-06-07 06:47:51 +08:00
|
|
|
|
2009-04-30 00:54:06 +08:00
|
|
|
if (masstotal > 0.0) return sqrt(rg_all/masstotal);
|
|
|
|
return 0.0;
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
2010-05-19 08:11:00 +08:00
|
|
|
compute the angular momentum L (lmom) of group
|
|
|
|
around center-of-mass cm
|
2006-09-28 03:51:33 +08:00
|
|
|
must unwrap atoms to compute L correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::angmom(int igroup, double *cm, double *lmom)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **v = atom->v;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2006-09-28 03:51:33 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz,massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
double p[3];
|
|
|
|
p[0] = p[1] = p[2] = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) massone = rmass[i];
|
|
|
|
else massone = mass[type[i]];
|
2006-09-28 03:51:33 +08:00
|
|
|
p[0] += massone * (dy*v[i][2] - dz*v[i][1]);
|
|
|
|
p[1] += massone * (dz*v[i][0] - dx*v[i][2]);
|
|
|
|
p[2] += massone * (dx*v[i][1] - dy*v[i][0]);
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(p,lmom,3,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
2010-05-19 08:11:00 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the angular momentum L (lmom) of group of atoms in region
|
|
|
|
around center-of-mass cm
|
|
|
|
must unwrap atoms to compute L correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::angmom(int igroup, double *cm, double *lmom, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2010-05-19 08:11:00 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **v = atom->v;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2010-05-19 08:11:00 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz,massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
|
|
|
|
2010-05-19 08:11:00 +08:00
|
|
|
double p[3];
|
|
|
|
p[0] = p[1] = p[2] = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2010-05-19 08:11:00 +08:00
|
|
|
if (rmass) massone = rmass[i];
|
|
|
|
else massone = mass[type[i]];
|
|
|
|
p[0] += massone * (dy*v[i][2] - dz*v[i][1]);
|
|
|
|
p[1] += massone * (dz*v[i][0] - dx*v[i][2]);
|
|
|
|
p[2] += massone * (dx*v[i][1] - dy*v[i][0]);
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(p,lmom,3,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
2010-11-24 22:39:50 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the torque T (tq) on group
|
|
|
|
around center-of-mass cm
|
|
|
|
must unwrap atoms to compute T correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::torque(int igroup, double *cm, double *tq)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **f = atom->f;
|
|
|
|
int *mask = atom->mask;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2010-11-24 22:39:50 +08:00
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
|
|
|
|
2010-11-24 22:39:50 +08:00
|
|
|
double tlocal[3];
|
|
|
|
tlocal[0] = tlocal[1] = tlocal[2] = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2010-11-24 22:39:50 +08:00
|
|
|
tlocal[0] += dy*f[i][2] - dz*f[i][1];
|
|
|
|
tlocal[1] += dz*f[i][0] - dx*f[i][2];
|
|
|
|
tlocal[2] += dx*f[i][1] - dy*f[i][0];
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(tlocal,tq,3,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute the torque T (tq) on group of atoms in region
|
|
|
|
around center-of-mass cm
|
|
|
|
must unwrap atoms to compute T correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::torque(int igroup, double *cm, double *tq, int iregion)
|
|
|
|
{
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2010-11-24 22:39:50 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
double **f = atom->f;
|
|
|
|
int *mask = atom->mask;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2010-11-24 22:39:50 +08:00
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
|
|
|
|
2010-11-24 22:39:50 +08:00
|
|
|
double tlocal[3];
|
|
|
|
tlocal[0] = tlocal[1] = tlocal[2] = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2010-11-24 22:39:50 +08:00
|
|
|
tlocal[0] += dy*f[i][2] - dz*f[i][1];
|
|
|
|
tlocal[1] += dz*f[i][0] - dx*f[i][2];
|
|
|
|
tlocal[2] += dx*f[i][1] - dy*f[i][0];
|
|
|
|
}
|
|
|
|
|
|
|
|
MPI_Allreduce(tlocal,tq,3,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
/* ----------------------------------------------------------------------
|
2010-05-29 00:44:34 +08:00
|
|
|
compute moment of inertia tensor around center-of-mass cm of group
|
2006-09-28 03:51:33 +08:00
|
|
|
must unwrap atoms to compute itensor correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::inertia(int igroup, double *cm, double itensor[3][3])
|
|
|
|
{
|
|
|
|
int i,j;
|
|
|
|
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2006-09-28 03:51:33 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz,massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
|
|
|
|
2006-09-28 03:51:33 +08:00
|
|
|
double ione[3][3];
|
|
|
|
for (i = 0; i < 3; i++)
|
|
|
|
for (j = 0; j < 3; j++)
|
|
|
|
ione[i][j] = 0.0;
|
|
|
|
|
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2009-02-13 05:36:52 +08:00
|
|
|
if (rmass) massone = rmass[i];
|
|
|
|
else massone = mass[type[i]];
|
2006-09-28 03:51:33 +08:00
|
|
|
ione[0][0] += massone * (dy*dy + dz*dz);
|
|
|
|
ione[1][1] += massone * (dx*dx + dz*dz);
|
|
|
|
ione[2][2] += massone * (dx*dx + dy*dy);
|
|
|
|
ione[0][1] -= massone * dx*dy;
|
2010-05-29 00:44:34 +08:00
|
|
|
ione[1][2] -= massone * dy*dz;
|
|
|
|
ione[0][2] -= massone * dx*dz;
|
|
|
|
}
|
|
|
|
ione[1][0] = ione[0][1];
|
|
|
|
ione[2][1] = ione[1][2];
|
|
|
|
ione[2][0] = ione[0][2];
|
|
|
|
|
|
|
|
MPI_Allreduce(&ione[0][0],&itensor[0][0],9,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
|
|
|
compute moment of inertia tensor around cm of group of atoms in region
|
|
|
|
must unwrap atoms to compute itensor correctly
|
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
|
|
|
void Group::inertia(int igroup, double *cm, double itensor[3][3], int iregion)
|
|
|
|
{
|
|
|
|
int i,j;
|
|
|
|
|
|
|
|
int groupbit = bitmask[igroup];
|
|
|
|
Region *region = domain->regions[iregion];
|
2014-05-02 22:28:56 +08:00
|
|
|
region->prematch();
|
2010-05-29 00:44:34 +08:00
|
|
|
|
|
|
|
double **x = atom->x;
|
|
|
|
int *mask = atom->mask;
|
|
|
|
int *type = atom->type;
|
2014-01-15 00:17:20 +08:00
|
|
|
imageint *image = atom->image;
|
2010-05-29 00:44:34 +08:00
|
|
|
double *mass = atom->mass;
|
|
|
|
double *rmass = atom->rmass;
|
|
|
|
int nlocal = atom->nlocal;
|
|
|
|
|
|
|
|
double dx,dy,dz,massone;
|
2012-10-26 00:24:37 +08:00
|
|
|
double unwrap[3];
|
|
|
|
|
2010-05-29 00:44:34 +08:00
|
|
|
double ione[3][3];
|
|
|
|
for (i = 0; i < 3; i++)
|
|
|
|
for (j = 0; j < 3; j++)
|
|
|
|
ione[i][j] = 0.0;
|
|
|
|
|
|
|
|
for (i = 0; i < nlocal; i++)
|
|
|
|
if (mask[i] & groupbit && region->match(x[i][0],x[i][1],x[i][2])) {
|
2012-10-26 00:24:37 +08:00
|
|
|
domain->unmap(x[i],image[i],unwrap);
|
|
|
|
dx = unwrap[0] - cm[0];
|
|
|
|
dy = unwrap[1] - cm[1];
|
|
|
|
dz = unwrap[2] - cm[2];
|
2010-05-29 00:44:34 +08:00
|
|
|
if (rmass) massone = rmass[i];
|
|
|
|
else massone = mass[type[i]];
|
|
|
|
ione[0][0] += massone * (dy*dy + dz*dz);
|
|
|
|
ione[1][1] += massone * (dx*dx + dz*dz);
|
|
|
|
ione[2][2] += massone * (dx*dx + dy*dy);
|
|
|
|
ione[0][1] -= massone * dx*dy;
|
2006-09-28 03:51:33 +08:00
|
|
|
ione[1][2] -= massone * dy*dz;
|
|
|
|
ione[0][2] -= massone * dx*dz;
|
|
|
|
}
|
|
|
|
ione[1][0] = ione[0][1];
|
|
|
|
ione[2][1] = ione[1][2];
|
|
|
|
ione[2][0] = ione[0][2];
|
|
|
|
|
|
|
|
MPI_Allreduce(&ione[0][0],&itensor[0][0],9,MPI_DOUBLE,MPI_SUM,world);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* ----------------------------------------------------------------------
|
2016-09-07 05:55:39 +08:00
|
|
|
compute angular velocity omega from L and I
|
2006-09-28 03:51:33 +08:00
|
|
|
------------------------------------------------------------------------- */
|
|
|
|
|
2007-01-30 08:22:05 +08:00
|
|
|
void Group::omega(double *angmom, double inertia[3][3], double *w)
|
2006-09-28 03:51:33 +08:00
|
|
|
{
|
2016-09-07 05:55:39 +08:00
|
|
|
double idiag[3],ex[3],ey[3],ez[3],cross[3];
|
2016-09-07 07:19:15 +08:00
|
|
|
double evectors[3][3],inverse[3][3];
|
|
|
|
|
|
|
|
// determinant = triple product of rows of inertia matrix
|
|
|
|
|
|
|
|
double determinant = inertia[0][0]*inertia[1][1]*inertia[2][2] +
|
|
|
|
inertia[0][1]*inertia[1][2]*inertia[2][0] +
|
|
|
|
inertia[0][2]*inertia[1][0]*inertia[2][1] -
|
|
|
|
inertia[0][0]*inertia[1][2]*inertia[2][1] -
|
|
|
|
inertia[0][1]*inertia[1][0]*inertia[2][2] -
|
|
|
|
inertia[2][0]*inertia[1][1]*inertia[0][2];
|
|
|
|
|
|
|
|
// non-singular I matrix
|
|
|
|
// use L = Iw, inverting I to solve for w
|
2016-09-08 00:12:51 +08:00
|
|
|
// this should give exact zeroing of angular momentum by velocity command
|
2016-09-07 07:19:15 +08:00
|
|
|
|
|
|
|
if (determinant > EPSILON) {
|
|
|
|
|
|
|
|
inverse[0][0] = inertia[1][1]*inertia[2][2] - inertia[1][2]*inertia[2][1];
|
|
|
|
inverse[0][1] = -(inertia[0][1]*inertia[2][2] -
|
|
|
|
inertia[0][2]*inertia[2][1]);
|
|
|
|
inverse[0][2] = inertia[0][1]*inertia[1][2] - inertia[0][2]*inertia[1][1];
|
|
|
|
|
|
|
|
inverse[1][0] = -(inertia[1][0]*inertia[2][2] -
|
|
|
|
inertia[1][2]*inertia[2][0]);
|
|
|
|
inverse[1][1] = inertia[0][0]*inertia[2][2] - inertia[0][2]*inertia[2][0];
|
|
|
|
inverse[1][2] = -(inertia[0][0]*inertia[1][2] -
|
|
|
|
inertia[0][2]*inertia[1][0]);
|
|
|
|
|
|
|
|
inverse[2][0] = inertia[1][0]*inertia[2][1] - inertia[1][1]*inertia[2][0];
|
|
|
|
inverse[2][1] = -(inertia[0][0]*inertia[2][1] -
|
|
|
|
inertia[0][1]*inertia[2][0]);
|
|
|
|
inverse[2][2] = inertia[0][0]*inertia[1][1] - inertia[0][1]*inertia[1][0];
|
|
|
|
|
|
|
|
for (int i = 0; i < 3; i++)
|
|
|
|
for (int j = 0; j < 3; j++)
|
|
|
|
inverse[i][j] /= determinant;
|
|
|
|
|
|
|
|
w[0] = inverse[0][0]*angmom[0] + inverse[0][1]*angmom[1] +
|
|
|
|
inverse[0][2]*angmom[2];
|
|
|
|
w[1] = inverse[1][0]*angmom[0] + inverse[1][1]*angmom[1] +
|
|
|
|
inverse[1][2]*angmom[2];
|
|
|
|
w[2] = inverse[2][0]*angmom[0] + inverse[2][1]*angmom[1] +
|
|
|
|
inverse[2][2]*angmom[2];
|
|
|
|
|
|
|
|
// handle (nearly) singular I matrix
|
2016-09-07 23:49:34 +08:00
|
|
|
// typically due to 2-atom group or linear molecule
|
2016-09-07 07:19:15 +08:00
|
|
|
// use jacobi() and angmom_to_omega() to calculate valid omega
|
2016-09-08 00:12:51 +08:00
|
|
|
// less exact answer than matrix inversion, due to iterative Jacobi method
|
2016-09-07 07:19:15 +08:00
|
|
|
|
|
|
|
} else {
|
|
|
|
int ierror = MathExtra::jacobi(inertia,idiag,evectors);
|
|
|
|
if (ierror) error->all(FLERR,
|
|
|
|
"Insufficient Jacobi rotations for group::omega");
|
|
|
|
|
|
|
|
ex[0] = evectors[0][0];
|
|
|
|
ex[1] = evectors[1][0];
|
|
|
|
ex[2] = evectors[2][0];
|
|
|
|
ey[0] = evectors[0][1];
|
|
|
|
ey[1] = evectors[1][1];
|
|
|
|
ey[2] = evectors[2][1];
|
|
|
|
ez[0] = evectors[0][2];
|
|
|
|
ez[1] = evectors[1][2];
|
|
|
|
ez[2] = evectors[2][2];
|
2016-09-07 05:55:39 +08:00
|
|
|
|
2016-09-07 07:19:15 +08:00
|
|
|
// enforce 3 evectors as a right-handed coordinate system
|
|
|
|
// flip 3rd vector if needed
|
2016-09-07 05:55:39 +08:00
|
|
|
|
2016-09-07 07:19:15 +08:00
|
|
|
MathExtra::cross3(ex,ey,cross);
|
|
|
|
if (MathExtra::dot3(cross,ez) < 0.0) MathExtra::negate3(ez);
|
2016-09-07 05:55:39 +08:00
|
|
|
|
2016-09-07 07:19:15 +08:00
|
|
|
// if any principal moment < scaled EPSILON, set to 0.0
|
2016-09-07 05:55:39 +08:00
|
|
|
|
2016-09-07 07:19:15 +08:00
|
|
|
double max;
|
|
|
|
max = MAX(idiag[0],idiag[1]);
|
|
|
|
max = MAX(max,idiag[2]);
|
|
|
|
|
|
|
|
if (idiag[0] < EPSILON*max) idiag[0] = 0.0;
|
|
|
|
if (idiag[1] < EPSILON*max) idiag[1] = 0.0;
|
|
|
|
if (idiag[2] < EPSILON*max) idiag[2] = 0.0;
|
|
|
|
|
|
|
|
// calculate omega using diagonalized inertia matrix
|
|
|
|
|
|
|
|
MathExtra::angmom_to_omega(angmom,ex,ey,ez,idiag,w);
|
|
|
|
}
|
2006-09-28 03:51:33 +08:00
|
|
|
}
|