modernize parsing of the Atoms section

This commit is contained in:
Axel Kohlmeyer
2021-12-29 20:24:27 -05:00
parent c97483c46f
commit 64d6a2fd1f
3 changed files with 19 additions and 26 deletions

View File

@ -1058,6 +1058,7 @@ void Atom::data_atoms(int n, char *buf, tagint id_offset, tagint mol_offset,
double *coord; double *coord;
char *next; char *next;
// use the first line to detect and validate the number of words/tokens per line
next = strchr(buf,'\n'); next = strchr(buf,'\n');
*next = '\0'; *next = '\0';
int nwords = utils::trim_and_count_words(buf); int nwords = utils::trim_and_count_words(buf);
@ -1066,8 +1067,6 @@ void Atom::data_atoms(int n, char *buf, tagint id_offset, tagint mol_offset,
if (nwords != avec->size_data_atom && nwords != avec->size_data_atom + 3) if (nwords != avec->size_data_atom && nwords != avec->size_data_atom + 3)
error->all(FLERR,"Incorrect atom format in data file"); error->all(FLERR,"Incorrect atom format in data file");
char **values = new char*[nwords];
// set bounds for my proc // set bounds for my proc
// if periodic and I am lo/hi proc, adjust bounds by EPSILON // if periodic and I am lo/hi proc, adjust bounds by EPSILON
// insures all data atoms will be owned even with round-off // insures all data atoms will be owned even with round-off
@ -1138,21 +1137,16 @@ void Atom::data_atoms(int n, char *buf, tagint id_offset, tagint mol_offset,
for (int i = 0; i < n; i++) { for (int i = 0; i < n; i++) {
next = strchr(buf,'\n'); next = strchr(buf,'\n');
*next = '\0';
for (m = 0; m < nwords; m++) { auto values = Tokenizer(utils::trim_comment(buf)).as_vector();
buf += strspn(buf," \t\n\r\f"); if (values.size() != nwords)
buf[strcspn(buf," \t\n\r\f")] = '\0'; error->all(FLERR, "Incorrect atom format in data file: {}", utils::trim(buf));
if (strlen(buf) == 0)
error->all(FLERR,"Incorrect atom format in data file");
values[m] = buf;
buf += strlen(buf)+1;
}
int imx = 0, imy = 0, imz = 0; int imx = 0, imy = 0, imz = 0;
if (imageflag) { if (imageflag) {
imx = utils::inumeric(FLERR,values[iptr],false,lmp); imx = utils::inumeric(FLERR,values[iptr].c_str(),false,lmp);
imy = utils::inumeric(FLERR,values[iptr+1],false,lmp); imy = utils::inumeric(FLERR,values[iptr+1].c_str(),false,lmp);
imz = utils::inumeric(FLERR,values[iptr+2],false,lmp); imz = utils::inumeric(FLERR,values[iptr+2].c_str(),false,lmp);
if ((domain->dimension == 2) && (imz != 0)) if ((domain->dimension == 2) && (imz != 0))
error->all(FLERR,"Z-direction image flag must be 0 for 2d-systems"); error->all(FLERR,"Z-direction image flag must be 0 for 2d-systems");
if ((!domain->xperiodic) && (imx != 0)) { reset_image_flag[0] = true; imx = 0; } if ((!domain->xperiodic) && (imx != 0)) { reset_image_flag[0] = true; imx = 0; }
@ -1163,9 +1157,9 @@ void Atom::data_atoms(int n, char *buf, tagint id_offset, tagint mol_offset,
(((imageint) (imy + IMGMAX) & IMGMASK) << IMGBITS) | (((imageint) (imy + IMGMAX) & IMGMASK) << IMGBITS) |
(((imageint) (imz + IMGMAX) & IMGMASK) << IMG2BITS); (((imageint) (imz + IMGMAX) & IMGMASK) << IMG2BITS);
xdata[0] = utils::numeric(FLERR,values[xptr],false,lmp); xdata[0] = utils::numeric(FLERR,values[xptr].c_str(),false,lmp);
xdata[1] = utils::numeric(FLERR,values[xptr+1],false,lmp); xdata[1] = utils::numeric(FLERR,values[xptr+1].c_str(),false,lmp);
xdata[2] = utils::numeric(FLERR,values[xptr+2],false,lmp); xdata[2] = utils::numeric(FLERR,values[xptr+2].c_str(),false,lmp);
if (shiftflag) { if (shiftflag) {
xdata[0] += shift[0]; xdata[0] += shift[0];
xdata[1] += shift[1]; xdata[1] += shift[1];
@ -1193,7 +1187,6 @@ void Atom::data_atoms(int n, char *buf, tagint id_offset, tagint mol_offset,
buf = next + 1; buf = next + 1;
} }
delete [] values;
} }
/* ---------------------------------------------------------------------- /* ----------------------------------------------------------------------

View File

@ -1707,7 +1707,7 @@ void AtomVec::create_atom(int itype, double *coord)
initialize other peratom quantities initialize other peratom quantities
------------------------------------------------------------------------- */ ------------------------------------------------------------------------- */
void AtomVec::data_atom(double *coord, imageint imagetmp, char **values) void AtomVec::data_atom(double *coord, imageint imagetmp, const std::vector<std::string> &values)
{ {
int m,n,datatype,cols; int m,n,datatype,cols;
void *pdata; void *pdata;
@ -1732,7 +1732,7 @@ void AtomVec::data_atom(double *coord, imageint imagetmp, char **values)
if (datatype == Atom::DOUBLE) { if (datatype == Atom::DOUBLE) {
if (cols == 0) { if (cols == 0) {
double *vec = *((double **) pdata); double *vec = *((double **) pdata);
vec[nlocal] = utils::numeric(FLERR,values[ivalue++],true,lmp); vec[nlocal] = utils::numeric(FLERR,values[ivalue++].c_str(),true,lmp);
} else { } else {
double **array = *((double ***) pdata); double **array = *((double ***) pdata);
if (array == atom->x) { // x was already set by coord arg if (array == atom->x) { // x was already set by coord arg
@ -1740,25 +1740,25 @@ void AtomVec::data_atom(double *coord, imageint imagetmp, char **values)
continue; continue;
} }
for (m = 0; m < cols; m++) for (m = 0; m < cols; m++)
array[nlocal][m] = utils::numeric(FLERR,values[ivalue++],true,lmp); array[nlocal][m] = utils::numeric(FLERR,values[ivalue++].c_str(),true,lmp);
} }
} else if (datatype == Atom::INT) { } else if (datatype == Atom::INT) {
if (cols == 0) { if (cols == 0) {
int *vec = *((int **) pdata); int *vec = *((int **) pdata);
vec[nlocal] = utils::inumeric(FLERR,values[ivalue++],true,lmp); vec[nlocal] = utils::inumeric(FLERR,values[ivalue++].c_str(),true,lmp);
} else { } else {
int **array = *((int ***) pdata); int **array = *((int ***) pdata);
for (m = 0; m < cols; m++) for (m = 0; m < cols; m++)
array[nlocal][m] = utils::inumeric(FLERR,values[ivalue++],true,lmp); array[nlocal][m] = utils::inumeric(FLERR,values[ivalue++].c_str(),true,lmp);
} }
} else if (datatype == Atom::BIGINT) { } else if (datatype == Atom::BIGINT) {
if (cols == 0) { if (cols == 0) {
bigint *vec = *((bigint **) pdata); bigint *vec = *((bigint **) pdata);
vec[nlocal] = utils::bnumeric(FLERR,values[ivalue++],true,lmp); vec[nlocal] = utils::bnumeric(FLERR,values[ivalue++].c_str(),true,lmp);
} else { } else {
bigint **array = *((bigint ***) pdata); bigint **array = *((bigint ***) pdata);
for (m = 0; m < cols; m++) for (m = 0; m < cols; m++)
array[nlocal][m] = utils::bnumeric(FLERR,values[ivalue++],true,lmp); array[nlocal][m] = utils::bnumeric(FLERR,values[ivalue++].c_str(),true,lmp);
} }
} }
} }

View File

@ -124,7 +124,7 @@ class AtomVec : protected Pointers {
virtual void create_atom(int, double *); virtual void create_atom(int, double *);
virtual void create_atom_post(int) {} virtual void create_atom_post(int) {}
virtual void data_atom(double *, imageint, char **); virtual void data_atom(double *, imageint, const std::vector<std::string> &);
virtual void data_atom_post(int) {} virtual void data_atom_post(int) {}
virtual void data_atom_bonus(int, char **) {} virtual void data_atom_bonus(int, char **) {}
virtual void data_body(int, int, int, int *, double *) {} virtual void data_body(int, int, int, int *, double *) {}