refactor group2ndx and ndx2group commands to use fmtlib, tokenizer and utils

This commit is contained in:
Axel Kohlmeyer
2021-03-28 13:12:39 -04:00
parent dfb18caf5a
commit 31726f56e6
6 changed files with 187 additions and 244 deletions

View File

@ -23,6 +23,8 @@
#include "error.h" #include "error.h"
#include "group.h" #include "group.h"
#include <cmath>
using namespace LAMMPS_NS; using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- /* ----------------------------------------------------------------------
@ -40,96 +42,6 @@ static int cmptagint(const void *p1, const void *p2)
} }
} }
/* ----------------------------------------------------------------------
helper function. writes out one group to a gromacs style index file
---------------------------------------------------------------------- */
static void write_group(FILE *fp, int gid, Atom *atom, Group *group, int me,
int np, MPI_Comm world, FILE *screen, FILE *logfile)
{
char fmt[16];
tagint *sendlist, *recvlist;
bigint num = group->count(gid);
int lnum, cols;
if (me == 0) {
if (screen) fprintf(screen, " writing group %s... ", group->names[gid]);
if (logfile) fprintf(logfile, " writing group %s... ", group->names[gid]);
// the "all" group in LAMMPS is called "System" in gromacs
if (gid == 0) {
fputs("[ System ]\n", fp);
} else {
fprintf(fp,"[ %s ]\n", group->names[gid]);
}
// derive format string for index lists
bigint j = atom->natoms;
int i=0;
while (j > 0) {
++i;
j /= 10;
}
snprintf(fmt,16,"%%%dd ", i);
cols = 80 / (i+1);
}
if (num > 0) {
const int * const mask = atom->mask;
const tagint * const tag = atom->tag;
const int groupbit = group->bitmask[gid];
const int nlocal = atom->nlocal;
int i;
sendlist = new tagint[nlocal];
recvlist = new tagint[num];
lnum = 0;
for (i = 0; i < nlocal; i++)
if (mask[i] & groupbit) sendlist[lnum++] = tag[i];
int nrecv,allrecv;
if (me == 0) {
MPI_Status status;
MPI_Request request;
for (i=0; i < lnum; i++)
recvlist[i] = sendlist[i];
allrecv = lnum;
for (i=1; i < np; ++i) {
MPI_Irecv(recvlist+allrecv,num-allrecv,MPI_LMP_TAGINT,i,0, world,&request);
MPI_Send(&nrecv,0,MPI_INT,i,0,world);
MPI_Wait(&request,&status);
MPI_Get_count(&status,MPI_LMP_TAGINT,&nrecv);
allrecv += nrecv;
}
// sort received list
qsort((void *)recvlist, num, sizeof(tagint), cmptagint);
} else {
MPI_Recv(&nrecv,0,MPI_INT,0,0,world,MPI_STATUS_IGNORE);
MPI_Rsend(sendlist,lnum,MPI_LMP_TAGINT,0,0,world);
}
delete [] sendlist;
}
if (me == 0) {
int i, j;
for (i=0, j=0; i < num; ++i) {
fprintf(fp,fmt,recvlist[i]);
++j;
if (j == cols) {
fputs("\n",fp);
j = 0;
}
}
if (j > 0) fputs("\n",fp);
if (screen) fputs("done\n",screen);
if (logfile) fputs("done\n",logfile);
}
if (num > 0) delete[] recvlist;
}
/* ---------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- */
void Group2Ndx::command(int narg, char **arg) void Group2Ndx::command(int narg, char **arg)
@ -144,31 +56,100 @@ void Group2Ndx::command(int narg, char **arg)
if (comm->me == 0) { if (comm->me == 0) {
fp = fopen(arg[0], "w"); fp = fopen(arg[0], "w");
if (fp == nullptr) if (fp == nullptr)
error->one(FLERR,"Cannot open index file for writing"); error->one(FLERR,fmt::format("Cannot open index file for writing: {}",
utils::getsyserror()));
if (screen) utils::logmesg(lmp,fmt::format("Writing groups to index file {}:\n",arg[0]));
fprintf(screen, "Writing groups to index file %s:\n",arg[0]);
if (logfile)
fprintf(logfile,"Writing groups to index file %s:\n",arg[0]);
} }
if (narg == 1) { // write out all groups if (narg == 1) { // write out all groups
for (int i=0; i < group->ngroup; ++i) { for (int i=0; i < group->ngroup; ++i) {
write_group(fp,i,atom,group,comm->me,comm->nprocs,world,screen,logfile); write_group(fp,i);
} }
} else { // write only selected groups } else { // write only selected groups
for (int i=1; i < narg; ++i) { for (int i=1; i < narg; ++i) {
int gid = group->find(arg[i]); int gid = group->find(arg[i]);
if (gid < 0) error->all(FLERR, "Non-existing group requested"); if (gid < 0) error->all(FLERR, "Non-existing group requested");
write_group(fp,gid,atom,group,comm->me,comm->nprocs,world,screen,logfile); write_group(fp,gid);
} }
} }
if (comm->me == 0) { if (comm->me == 0) fclose(fp);
if (screen) fputs("\n",screen);
if (logfile) fputs("\n",logfile);
fclose(fp);
}
} }
/* ----------------------------------------------------------------------
write out one group to a Gromacs style index file
---------------------------------------------------------------------- */
void Group2Ndx::write_group(FILE *fp, int gid)
{
tagint *sendlist, *recvlist;
bigint gcount = group->count(gid);
int lnum, width, cols;
if (comm->me == 0) {
utils::logmesg(lmp,fmt::format(" writing group {}...",group->names[gid]));
// the "all" group in LAMMPS is called "System" in Gromacs
if (gid == 0) {
fputs("[ System ]\n", fp);
} else {
fmt::print(fp,"[ {} ]\n", group->names[gid]);
}
width = log10((double) atom->natoms)+2;
cols = 80 / width;
}
if (gcount > 0) {
const int * const mask = atom->mask;
const tagint * const tag = atom->tag;
const int groupbit = group->bitmask[gid];
const int nlocal = atom->nlocal;
int i;
sendlist = new tagint[nlocal];
recvlist = new tagint[gcount];
lnum = 0;
for (i = 0; i < nlocal; i++)
if (mask[i] & groupbit) sendlist[lnum++] = tag[i];
int nrecv=0;
bigint allrecv;
if (comm->me == 0) {
MPI_Status status;
MPI_Request request;
for (i=0; i < lnum; i++)
recvlist[i] = sendlist[i];
allrecv = lnum;
for (i=1; i < comm->nprocs; ++i) {
MPI_Irecv(recvlist+allrecv,gcount-allrecv,MPI_LMP_TAGINT,i,0, world,&request);
MPI_Send(&nrecv,0,MPI_INT,i,0,world); // block rank "i" until we are ready to receive
MPI_Wait(&request,&status);
MPI_Get_count(&status,MPI_LMP_TAGINT,&nrecv);
allrecv += nrecv;
}
// sort received list
qsort((void *)recvlist, allrecv, sizeof(tagint), cmptagint);
} else {
MPI_Recv(&nrecv,0,MPI_INT,0,0,world,MPI_STATUS_IGNORE);
MPI_Rsend(sendlist,lnum,MPI_LMP_TAGINT,0,0,world);
}
delete [] sendlist;
}
if (comm->me == 0) {
int i, j;
for (i=0, j=0; i < gcount; ++i) {
fmt::print(fp,"{:>{}}",recvlist[i],width);
++j;
if (j == cols) {
fputs("\n",fp);
j = 0;
}
}
if (j > 0) fputs("\n",fp);
utils::logmesg(lmp,"done\n");
}
if (gcount > 0) delete[] recvlist;
}

View File

@ -30,6 +30,8 @@ class Group2Ndx : protected Pointers {
public: public:
Group2Ndx(class LAMMPS *lmp) : Pointers(lmp) {}; Group2Ndx(class LAMMPS *lmp) : Pointers(lmp) {};
void command(int, char **); void command(int, char **);
private:
void write_group(FILE *, int);
}; };
} }

View File

@ -22,64 +22,52 @@
#include "comm.h" #include "comm.h"
#include "error.h" #include "error.h"
#include "group.h" #include "group.h"
#include "memory.h"
#include "tokenizer.h"
#include <cstring> #include <cstdlib>
using namespace LAMMPS_NS; using namespace LAMMPS_NS;
#define BUFLEN 4096 #define BUFLEN 4096
#define DELTA 16384 #define DELTA 16384
static char *find_section(FILE *fp, const char *name) // read file until next section "name" or any next section if name == ""
static std::string find_section(FILE *fp, const std::string &name)
{ {
char linebuf[BUFLEN]; char linebuf[BUFLEN];
char *n,*p,*t,*r;
while ((p = fgets(linebuf,BUFLEN,fp))) { std::string pattern = "^\\s*\\[\\s+\\S+\\s+\\]\\s*$";
t = strtok(p," \t\n\r\f"); if (!name.empty())
if ((t != nullptr) && *t == '[') { pattern = fmt::format("^\\s*\\[\\s+{}\\s+\\]\\s*$",name);
t = strtok(nullptr," \t\n\r\f");
if (t != nullptr) { fgets(linebuf,BUFLEN,fp);
n = t; while (!feof(fp)) {
t = strtok(nullptr," \t\n\r\f"); if (utils::strmatch(linebuf,pattern))
if ((t != nullptr) && *t == ']') { return Tokenizer(linebuf).as_vector()[1];
if ((name == nullptr) || strcmp(name,n) == 0) { fgets(linebuf,BUFLEN,fp);
int l = strlen(n);
r = new char[l+1];
strncpy(r,n,l+1);
return r;
}
}
}
}
} }
return nullptr; return "";
} }
static tagint *read_section(FILE *fp, bigint &num) static std::vector<tagint> read_section(FILE *fp, std::string &name)
{ {
char linebuf[BUFLEN]; char linebuf[BUFLEN];
char *p,*t; std::vector<tagint> tagbuf;
tagint *tagbuf; std::string pattern = "^\\s*\\[\\s+\\S+\\s+\\]\\s*$";
bigint nmax;
num = 0; while (fgets(linebuf,BUFLEN,fp)) {
nmax = DELTA; // start of new section. we are done, update "name"
tagbuf = (tagint *)malloc(sizeof(tagint)*nmax); if (utils::strmatch(linebuf,pattern)) {
name = Tokenizer(linebuf).as_vector()[1];
while ((p = fgets(linebuf,BUFLEN,fp))) { return tagbuf;
t = strtok(p," \t\n\r\f");
while (t != nullptr) {
// start of a new section. we are done here.
if (*t == '[') return tagbuf;
tagbuf[num++] = ATOTAGINT(t);
if (num == nmax) {
nmax += DELTA;
tagbuf = (tagint *)realloc(tagbuf,sizeof(tagint)*nmax);
}
t = strtok(nullptr," \t\n\r\f");
} }
ValueTokenizer values(linebuf);
while (values.has_next())
tagbuf.push_back(values.next_tagint());
} }
// set empty name to indicate end of file
name = "";
return tagbuf; return tagbuf;
} }
@ -90,151 +78,122 @@ void Ndx2Group::command(int narg, char **arg)
int len; int len;
bigint num; bigint num;
FILE *fp; FILE *fp;
char *name = nullptr; tagint *tagbuf;
tagint *tags; std::string name = "", next;
if (narg < 1) error->all(FLERR,"Illegal ndx2group command"); if (narg < 1) error->all(FLERR,"Illegal ndx2group command");
if (atom->tag_enable == 0) if (atom->tag_enable == 0)
error->all(FLERR,"Must have atom IDs for ndx2group command"); error->all(FLERR,"Must have atom IDs for ndx2group command");
if (comm->me == 0) { if (comm->me == 0) {
fp = fopen(arg[0], "r"); fp = fopen(arg[0], "r");
if (fp == nullptr) if (fp == nullptr)
error->one(FLERR,"Cannot open index file for reading"); error->one(FLERR,fmt::format("Cannot open index file for reading: {}",
utils::getsyserror()));
if (screen) utils::logmesg(lmp,fmt::format("Reading groups from index file {}:\n",arg[0]));
fprintf(screen, "Reading groups from index file %s:\n",arg[0]);
if (logfile)
fprintf(logfile,"Reading groups from index file %s:\n",arg[0]);
} }
if (narg == 1) { // restore all groups if (narg == 1) { // restore all groups
do { if (comm->me == 0) {
if (comm->me == 0) { name = find_section(fp,"");
len = 0; while (!name.empty()) {
// find the next section. // skip over group "all", which is called "System" in gromacs
// if we had processed a section, before we need to step back if (name == "System") {
if (name != nullptr) { name = find_section(fp,"");
rewind(fp); continue;
char *tmp = find_section(fp,name);
delete[] tmp;
delete[] name;
name = nullptr;
} }
name = find_section(fp,nullptr);
if (name != nullptr) {
len=strlen(name)+1;
// skip over group "all", which is called "System" in gromacs utils::logmesg(lmp,fmt::format(" Processing group '{}'\n",name));
if (strcmp(name,"System") == 0) continue; len = name.size()+1;
if (screen)
fprintf(screen," Processing group '%s'\n",name);
if (logfile)
fprintf(logfile," Processing group '%s'\n",name);
}
MPI_Bcast(&len,1,MPI_INT,0,world); MPI_Bcast(&len,1,MPI_INT,0,world);
if (len > 0) { if (len > 1) {
MPI_Bcast(name,len,MPI_CHAR,0,world); MPI_Bcast((void *)name.c_str(),len,MPI_CHAR,0,world);
// read tags for atoms in group and broadcast // read tags for atoms in group and broadcast
num = 0; std::vector<tagint> tags = read_section(fp,next);
tags = read_section(fp,num); num = tags.size();
MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world); MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world);
MPI_Bcast(tags,num,MPI_LMP_TAGINT,0,world); MPI_Bcast((void *)tags.data(),num,MPI_LMP_TAGINT,0,world);
create(name,num,tags); create(name,tags);
free(tags); name = next;
}
} else {
MPI_Bcast(&len,1,MPI_INT,0,world);
if (len > 0) {
delete[] name;
name = new char[len];
MPI_Bcast(name,len,MPI_CHAR,0,world);
MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world);
tags = (tagint *)malloc(sizeof(tagint)*(num ? num : 1));
MPI_Bcast(tags,num,MPI_LMP_TAGINT,0,world);
create(name,num,tags);
free(tags);
} }
} }
} while (len); len = -1;
MPI_Bcast(&len,1,MPI_INT,0,world);
} else {
while (1) {
MPI_Bcast(&len,1,MPI_INT,0,world);
if (len < 0) break;
if (len > 1) {
char *buf = new char[len];
MPI_Bcast(buf,len,MPI_CHAR,0,world);
MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world);
tagint *tbuf = new tagint[num];
MPI_Bcast(tbuf,num,MPI_LMP_TAGINT,0,world);
create(buf,std::vector<tagint>(tbuf,tbuf+num));
delete[] buf;
delete[] tbuf;
}
}
}
} else { // restore selected groups } else { // restore selected groups
for (int idx=1; idx < narg; ++idx) {
for (int idx=1; idx < narg; ++idx) {
if (comm->me == 0) { if (comm->me == 0) {
len = 0;
// find named section, search from beginning of file // find named section, search from beginning of file
if (name != nullptr) delete[] name;
rewind(fp); rewind(fp);
name = find_section(fp,arg[idx]); name = find_section(fp,arg[idx]);
if (name != nullptr) len=strlen(name)+1; utils::logmesg(lmp,fmt::format(" {} group '{}'\n", name.size()
? "Processing" : "Skipping",arg[idx]));
if (screen) len = name.size()+1;
fprintf(screen," %s group '%s'\n",
len ? "Processing" : "Skipping",arg[idx]);
if (logfile)
fprintf(logfile,"%s group '%s'\n",
len ? "Processing" : "Skipping",arg[idx]);
MPI_Bcast(&len,1,MPI_INT,0,world); MPI_Bcast(&len,1,MPI_INT,0,world);
if (len > 0) { if (len > 1) {
MPI_Bcast(name,len,MPI_CHAR,0,world); MPI_Bcast((void *)name.c_str(),len,MPI_CHAR,0,world);
// read tags for atoms in group and broadcast // read tags for atoms in group and broadcast
num = 0; std::vector<tagint> tags = read_section(fp,name);
tags = read_section(fp,num); num = tags.size();
MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world); MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world);
MPI_Bcast(tags,num,MPI_LMP_TAGINT,0,world); MPI_Bcast((void *)tags.data(),num,MPI_LMP_TAGINT,0,world);
create(name,num,tags); create(name,tags);
free(tags);
} }
} else { } else {
MPI_Bcast(&len,1,MPI_INT,0,world); MPI_Bcast(&len,1,MPI_INT,0,world);
if (len > 0) { if (len > 1) {
delete[] name; char *buf = new char[len];
name = new char[len]; MPI_Bcast(buf,len,MPI_CHAR,0,world);
MPI_Bcast(name,len,MPI_CHAR,0,world);
MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world); MPI_Bcast(&num,1,MPI_LMP_BIGINT,0,world);
tags = (tagint *)malloc(sizeof(tagint)*(num ? num : 1)); tagint *tbuf = new tagint[num];
MPI_Bcast(tags,num,MPI_LMP_TAGINT,0,world); MPI_Bcast(tbuf,num,MPI_LMP_TAGINT,0,world);
create(name,num,tags); create(buf,std::vector<tagint>(tbuf,tbuf+num));
free(tags); delete[] buf;
delete[] tbuf;
} }
} }
} }
} }
if (comm->me == 0) fclose(fp);
delete[] name;
if (comm->me == 0) {
if (screen) fputs("\n",screen);
if (logfile) fputs("\n",logfile);
fclose(fp);
}
} }
/* ---------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- */
void Ndx2Group::create(char *name, bigint num, tagint *tags) void Ndx2Group::create(const std::string &name, const std::vector<tagint> &tags)
{ {
// wipe out all members if the group exists. gid==0 is group "all" // wipe out all members if the group exists. gid==0 is group "all"
int gid = group->find(name); int gid = group->find(name);
if (gid > 0) group->assign(std::string(name) + " clear"); if (gid > 0) group->assign(name + " clear");
// map from global to local // map from global to local
const int nlocal = atom->nlocal; const int nlocal = atom->nlocal;
int *flags = (int *)calloc(nlocal,sizeof(int)); int *flags = (int *)calloc(nlocal,sizeof(int));
for (bigint i=0; i < num; ++i) { for (bigint i=0; i < tags.size(); ++i) {
const int id = atom->map(tags[i]); const int id = atom->map(tags[i]);
if (id < nlocal && id >= 0) if (id < nlocal && id >= 0) flags[id] = 1;
flags[id] = 1;
} }
group->create(name,flags); group->create(name,flags);
free(flags); free(flags);

View File

@ -23,6 +23,7 @@ CommandStyle(ndx2group,Ndx2Group)
#define LMP_NDX_GROUP_H #define LMP_NDX_GROUP_H
#include "pointers.h" #include "pointers.h"
#include <vector>
namespace LAMMPS_NS { namespace LAMMPS_NS {
@ -30,7 +31,7 @@ class Ndx2Group : protected Pointers {
public: public:
Ndx2Group(class LAMMPS *lmp) : Pointers(lmp) {}; Ndx2Group(class LAMMPS *lmp) : Pointers(lmp) {};
void command(int, char **); void command(int, char **);
void create(char *, bigint, tagint *); void create(const std::string &, const std::vector<tagint> &);
}; };
} }

View File

@ -553,7 +553,7 @@ void Group::assign(const std::string &groupcmd)
add flagged atoms to a new or existing group add flagged atoms to a new or existing group
------------------------------------------------------------------------- */ ------------------------------------------------------------------------- */
void Group::create(const char *name, int *flag) void Group::create(const std::string &name, int *flag)
{ {
int i; int i;

View File

@ -32,7 +32,7 @@ class Group : protected Pointers {
~Group(); ~Group();
void assign(int, char **); // assign atoms to a group void assign(int, char **); // assign atoms to a group
void assign(const std::string &); // convenience function void assign(const std::string &); // convenience function
void create(const char *, int *); // add flagged atoms to a group void create(const std::string &, int *); // add flagged atoms to a group
int find(const std::string &); // lookup name in list of groups int find(const std::string &); // lookup name in list of groups
int find_or_create(const char *); // lookup name or create new group int find_or_create(const char *); // lookup name or create new group
void write_restart(FILE *); void write_restart(FILE *);