add mechanism to check for known data file section names

using this mechanism we can reject custom section names that will
conflict with existing section names and thus avoid misleading errors.
apply this also to fix property atom, where the section name is
determined by the fix ID.
in addition, allow to specify NULL as section name, which will use
the fix ID.
This commit is contained in:
Axel Kohlmeyer
2021-10-28 14:23:27 -04:00
parent 4a048e3f57
commit adf1beea74
4 changed files with 76 additions and 68 deletions

View File

@ -254,7 +254,9 @@ the fix defines the syntax of the header line(s) and section that it
reads from the data file. Note that the *header-string* can be reads from the data file. Note that the *header-string* can be
specified as NULL, in which case no header lines are passed to the specified as NULL, in which case no header lines are passed to the
fix. This means the fix can infer the length of its Section from fix. This means the fix can infer the length of its Section from
standard header settings, such as the number of atoms. standard header settings, such as the number of atoms. Also the
*section-string* may be specified as NULL, and in that case the fix
ID is used as section name.
The formatting of individual lines in the data file (indentation, The formatting of individual lines in the data file (indentation,
spacing between words and numbers) is not important except that header spacing between words and numbers) is not important except that header

View File

@ -18,6 +18,7 @@
#include "comm.h" #include "comm.h"
#include "error.h" #include "error.h"
#include "memory.h" #include "memory.h"
#include "read_data.h"
#include "tokenizer.h" #include "tokenizer.h"
#include <cstring> #include <cstring>
@ -54,8 +55,7 @@ FixPropertyAtom::FixPropertyAtom(LAMMPS *lmp, int narg, char **arg) :
while (iarg < narg) { while (iarg < narg) {
if (strcmp(arg[iarg],"mol") == 0) { if (strcmp(arg[iarg],"mol") == 0) {
if (atom->molecule_flag) if (atom->molecule_flag)
error->all(FLERR,"Fix property/atom mol when atom_style " error->all(FLERR,"Fix property/atom mol when atom_style already has molecule attribute");
"already has molecule attribute");
if (molecule_flag) if (molecule_flag)
error->all(FLERR,"Fix property/atom cannot specify mol twice"); error->all(FLERR,"Fix property/atom cannot specify mol twice");
styles[nvalue] = MOLECULE; styles[nvalue] = MOLECULE;
@ -95,6 +95,8 @@ FixPropertyAtom::FixPropertyAtom(LAMMPS *lmp, int narg, char **arg) :
index[nvalue] = atom->find_custom(&arg[iarg][2],flag,ncols); index[nvalue] = atom->find_custom(&arg[iarg][2],flag,ncols);
if (index[nvalue] >= 0) if (index[nvalue] >= 0)
error->all(FLERR,"Fix property/atom vector name already exists"); error->all(FLERR,"Fix property/atom vector name already exists");
if (ReadData::is_data_section(id))
error->all(FLERR,"Fix property/atom fix ID must not be a data file section name");
index[nvalue] = atom->add_custom(&arg[iarg][2],0,0); index[nvalue] = atom->add_custom(&arg[iarg][2],0,0);
cols[nvalue] = 0; cols[nvalue] = 0;
values_peratom++; values_peratom++;
@ -107,6 +109,8 @@ FixPropertyAtom::FixPropertyAtom(LAMMPS *lmp, int narg, char **arg) :
index[nvalue] = atom->find_custom(&arg[iarg][2],flag,ncols); index[nvalue] = atom->find_custom(&arg[iarg][2],flag,ncols);
if (index[nvalue] >= 0) if (index[nvalue] >= 0)
error->all(FLERR,"Fix property/atom vector name already exists"); error->all(FLERR,"Fix property/atom vector name already exists");
if (ReadData::is_data_section(id))
error->all(FLERR,"Fix property/atom fix ID must not be a data file section name");
index[nvalue] = atom->add_custom(&arg[iarg][2],1,0); index[nvalue] = atom->add_custom(&arg[iarg][2],1,0);
cols[nvalue] = 0; cols[nvalue] = 0;
values_peratom++; values_peratom++;
@ -122,6 +126,8 @@ FixPropertyAtom::FixPropertyAtom(LAMMPS *lmp, int narg, char **arg) :
which = atom->find_custom(&arg[iarg][3],flag,ncols); which = atom->find_custom(&arg[iarg][3],flag,ncols);
if (which >= 0) if (which >= 0)
error->all(FLERR,"Fix property/atom array name {} already exists", &arg[iarg][3]); error->all(FLERR,"Fix property/atom array name {} already exists", &arg[iarg][3]);
if (ReadData::is_data_section(id))
error->all(FLERR,"Fix property/atom fix ID must not be a data file section name");
ncols = utils::inumeric(FLERR,arg[iarg+1],true,lmp); ncols = utils::inumeric(FLERR,arg[iarg+1],true,lmp);
if (ncols < 1) if (ncols < 1)

View File

@ -39,6 +39,8 @@
#include <cctype> #include <cctype>
#include <cstring> #include <cstring>
#include <string>
#include <unordered_set>
using namespace LAMMPS_NS; using namespace LAMMPS_NS;
@ -49,9 +51,27 @@ static constexpr int DELTA = 4; // must be 2 or larger
static constexpr int MAXBODY = 32; // max # of lines in one body static constexpr int MAXBODY = 32; // max # of lines in one body
// customize for new sections // customize for new sections
// change when add to header::section_keywords
static constexpr int NSECTIONS = 25;
static std::unordered_set<std::string> section_keywords = {
"Atoms", "Velocities", "Ellipsoids", "Lines", "Triangles", "Bodies",
"Bonds", "Angles", "Dihedrals", "Impropers",
"Masses", "Pair Coeffs", "PairIJ Coeffs", "Bond Coeffs", "Angle Coeffs",
"Dihedral Coeffs", "Improper Coeffs",
"BondBond Coeffs", "BondAngle Coeffs", "MiddleBondTorsion Coeffs",
"EndBondTorsion Coeffs", "AngleTorsion Coeffs",
"AngleAngleTorsion Coeffs", "BondBond13 Coeffs", "AngleAngle Coeffs"
};
// function to check whether a string is a known data section name
// made a static class member, so it can be called from other classes
bool ReadData::is_data_section(const std::string &keyword)
{
return section_keywords.count(keyword) > 0;
}
// clang-format off
enum{NONE,APPEND,VALUE,MERGE}; enum{NONE,APPEND,VALUE,MERGE};
// pair style suffixes to ignore // pair style suffixes to ignore
@ -266,18 +286,19 @@ void ReadData::command(int narg, char **arg)
if (iarg+4 > narg) if (iarg+4 > narg)
error->all(FLERR,"Illegal read_data command"); error->all(FLERR,"Illegal read_data command");
memory->grow(fix_index,nfix+1,"read_data:fix_index"); memory->grow(fix_index,nfix+1,"read_data:fix_index");
fix_header = (char **) fix_header = (char **) memory->srealloc(fix_header,(nfix+1)*sizeof(char *),
memory->srealloc(fix_header,(nfix+1)*sizeof(char *), "read_data:fix_header");
"read_data:fix_header"); fix_section = (char **) memory->srealloc(fix_section,(nfix+1)*sizeof(char *),
fix_section = (char **) "read_data:fix_section");
memory->srealloc(fix_section,(nfix+1)*sizeof(char *), if (is_data_section(arg[iarg+3]))
"read_data:fix_section"); error->all(FLERR,"Custom data section name {} for fix {} collides with existing "
"data section",arg[iarg+3],arg[iarg+1]);
fix_index[nfix] = modify->find_fix(arg[iarg+1]); fix_index[nfix] = modify->find_fix(arg[iarg+1]);
if (fix_index[nfix] < 0) if (fix_index[nfix] < 0) error->all(FLERR,"Fix ID for read_data does not exist");
error->all(FLERR,"Fix ID for read_data does not exist");
if (strcmp(arg[iarg+2],"NULL") == 0) fix_header[nfix] = nullptr; if (strcmp(arg[iarg+2],"NULL") == 0) fix_header[nfix] = nullptr;
else fix_header[nfix] = utils::strdup(arg[iarg+2]); else fix_header[nfix] = utils::strdup(arg[iarg+2]);
fix_section[nfix] = utils::strdup(arg[iarg+3]); if (strcmp(arg[iarg+3],"NULL") == 0) fix_section[nfix] = utils::strdup(arg[iarg+1]);
else fix_section[nfix] = utils::strdup(arg[iarg+3]);
nfix++; nfix++;
iarg += 4; iarg += 4;
@ -521,6 +542,7 @@ void ReadData::command(int narg, char **arg)
"from currently defined atom style"); "from currently defined atom style");
atoms(); atoms();
} else skip_lines(natoms); } else skip_lines(natoms);
} else if (strcmp(keyword,"Velocities") == 0) { } else if (strcmp(keyword,"Velocities") == 0) {
if (atomflag == 0) if (atomflag == 0)
error->all(FLERR,"Must read Atoms before Velocities"); error->all(FLERR,"Must read Atoms before Velocities");
@ -529,52 +551,50 @@ void ReadData::command(int narg, char **arg)
} else if (strcmp(keyword,"Bonds") == 0) { } else if (strcmp(keyword,"Bonds") == 0) {
topoflag = bondflag = 1; topoflag = bondflag = 1;
if (nbonds == 0) if (nbonds == 0) error->all(FLERR,"Invalid data file section: Bonds");
error->all(FLERR,"Invalid data file section: Bonds");
if (atomflag == 0) error->all(FLERR,"Must read Atoms before Bonds"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Bonds");
bonds(firstpass); bonds(firstpass);
} else if (strcmp(keyword,"Angles") == 0) { } else if (strcmp(keyword,"Angles") == 0) {
topoflag = angleflag = 1; topoflag = angleflag = 1;
if (nangles == 0) if (nangles == 0) error->all(FLERR,"Invalid data file section: Angles");
error->all(FLERR,"Invalid data file section: Angles");
if (atomflag == 0) error->all(FLERR,"Must read Atoms before Angles"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Angles");
angles(firstpass); angles(firstpass);
} else if (strcmp(keyword,"Dihedrals") == 0) { } else if (strcmp(keyword,"Dihedrals") == 0) {
topoflag = dihedralflag = 1; topoflag = dihedralflag = 1;
if (ndihedrals == 0) if (ndihedrals == 0) error->all(FLERR,"Invalid data file section: Dihedrals");
error->all(FLERR,"Invalid data file section: Dihedrals");
if (atomflag == 0) error->all(FLERR,"Must read Atoms before Dihedrals"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Dihedrals");
dihedrals(firstpass); dihedrals(firstpass);
} else if (strcmp(keyword,"Impropers") == 0) { } else if (strcmp(keyword,"Impropers") == 0) {
topoflag = improperflag = 1; topoflag = improperflag = 1;
if (nimpropers == 0) if (nimpropers == 0) error->all(FLERR,"Invalid data file section: Impropers");
error->all(FLERR,"Invalid data file section: Impropers");
if (atomflag == 0) error->all(FLERR,"Must read Atoms before Impropers"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Impropers");
impropers(firstpass); impropers(firstpass);
} else if (strcmp(keyword,"Ellipsoids") == 0) { } else if (strcmp(keyword,"Ellipsoids") == 0) {
ellipsoidflag = 1; ellipsoidflag = 1;
if (!avec_ellipsoid) if (!avec_ellipsoid) error->all(FLERR,"Invalid data file section: Ellipsoids");
error->all(FLERR,"Invalid data file section: Ellipsoids"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Ellipsoids");
if (atomflag == 0)
error->all(FLERR,"Must read Atoms before Ellipsoids");
if (firstpass) if (firstpass)
bonus(nellipsoids,(AtomVec *) avec_ellipsoid,"ellipsoids"); bonus(nellipsoids,(AtomVec *) avec_ellipsoid,"ellipsoids");
else skip_lines(nellipsoids); else skip_lines(nellipsoids);
} else if (strcmp(keyword,"Lines") == 0) { } else if (strcmp(keyword,"Lines") == 0) {
lineflag = 1; lineflag = 1;
if (!avec_line) if (!avec_line) error->all(FLERR,"Invalid data file section: Lines");
error->all(FLERR,"Invalid data file section: Lines");
if (atomflag == 0) error->all(FLERR,"Must read Atoms before Lines"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Lines");
if (firstpass) bonus(nlines,(AtomVec *) avec_line,"lines"); if (firstpass) bonus(nlines,(AtomVec *) avec_line,"lines");
else skip_lines(nlines); else skip_lines(nlines);
} else if (strcmp(keyword,"Triangles") == 0) { } else if (strcmp(keyword,"Triangles") == 0) {
triflag = 1; triflag = 1;
if (!avec_tri) if (!avec_tri) error->all(FLERR,"Invalid data file section: Triangles");
error->all(FLERR,"Invalid data file section: Triangles");
if (atomflag == 0) error->all(FLERR,"Must read Atoms before Triangles"); if (atomflag == 0) error->all(FLERR,"Must read Atoms before Triangles");
if (firstpass) bonus(ntris,(AtomVec *) avec_tri,"triangles"); if (firstpass) bonus(ntris,(AtomVec *) avec_tri,"triangles");
else skip_lines(ntris); else skip_lines(ntris);
} else if (strcmp(keyword,"Bodies") == 0) { } else if (strcmp(keyword,"Bodies") == 0) {
bodyflag = 1; bodyflag = 1;
if (!avec_body) if (!avec_body)
@ -655,6 +675,7 @@ void ReadData::command(int narg, char **arg)
error->all(FLERR,"Must define angle_style before BondBond Coeffs"); error->all(FLERR,"Must define angle_style before BondBond Coeffs");
if (firstpass) anglecoeffs(1); if (firstpass) anglecoeffs(1);
else skip_lines(nangletypes); else skip_lines(nangletypes);
} else if (strcmp(keyword,"BondAngle Coeffs") == 0) { } else if (strcmp(keyword,"BondAngle Coeffs") == 0) {
if (atom->avec->angles_allow == 0) if (atom->avec->angles_allow == 0)
error->all(FLERR,"Invalid data file section: BondAngle Coeffs"); error->all(FLERR,"Invalid data file section: BondAngle Coeffs");
@ -665,46 +686,41 @@ void ReadData::command(int narg, char **arg)
} else if (strcmp(keyword,"MiddleBondTorsion Coeffs") == 0) { } else if (strcmp(keyword,"MiddleBondTorsion Coeffs") == 0) {
if (atom->avec->dihedrals_allow == 0) if (atom->avec->dihedrals_allow == 0)
error->all(FLERR, error->all(FLERR,"Invalid data file section: MiddleBondTorsion Coeffs");
"Invalid data file section: MiddleBondTorsion Coeffs");
if (force->dihedral == nullptr) if (force->dihedral == nullptr)
error->all(FLERR, error->all(FLERR,"Must define dihedral_style before MiddleBondTorsion Coeffs");
"Must define dihedral_style before "
"MiddleBondTorsion Coeffs");
if (firstpass) dihedralcoeffs(1); if (firstpass) dihedralcoeffs(1);
else skip_lines(ndihedraltypes); else skip_lines(ndihedraltypes);
} else if (strcmp(keyword,"EndBondTorsion Coeffs") == 0) { } else if (strcmp(keyword,"EndBondTorsion Coeffs") == 0) {
if (atom->avec->dihedrals_allow == 0) if (atom->avec->dihedrals_allow == 0)
error->all(FLERR,"Invalid data file section: EndBondTorsion Coeffs"); error->all(FLERR,"Invalid data file section: EndBondTorsion Coeffs");
if (force->dihedral == nullptr) if (force->dihedral == nullptr)
error->all(FLERR, error->all(FLERR,"Must define dihedral_style before EndBondTorsion Coeffs");
"Must define dihedral_style before EndBondTorsion Coeffs");
if (firstpass) dihedralcoeffs(2); if (firstpass) dihedralcoeffs(2);
else skip_lines(ndihedraltypes); else skip_lines(ndihedraltypes);
} else if (strcmp(keyword,"AngleTorsion Coeffs") == 0) { } else if (strcmp(keyword,"AngleTorsion Coeffs") == 0) {
if (atom->avec->dihedrals_allow == 0) if (atom->avec->dihedrals_allow == 0)
error->all(FLERR,"Invalid data file section: AngleTorsion Coeffs"); error->all(FLERR,"Invalid data file section: AngleTorsion Coeffs");
if (force->dihedral == nullptr) if (force->dihedral == nullptr)
error->all(FLERR, error->all(FLERR,"Must define dihedral_style before AngleTorsion Coeffs");
"Must define dihedral_style before AngleTorsion Coeffs");
if (firstpass) dihedralcoeffs(3); if (firstpass) dihedralcoeffs(3);
else skip_lines(ndihedraltypes); else skip_lines(ndihedraltypes);
} else if (strcmp(keyword,"AngleAngleTorsion Coeffs") == 0) { } else if (strcmp(keyword,"AngleAngleTorsion Coeffs") == 0) {
if (atom->avec->dihedrals_allow == 0) if (atom->avec->dihedrals_allow == 0)
error->all(FLERR, error->all(FLERR,"Invalid data file section: AngleAngleTorsion Coeffs");
"Invalid data file section: AngleAngleTorsion Coeffs");
if (force->dihedral == nullptr) if (force->dihedral == nullptr)
error->all(FLERR, error->all(FLERR,"Must define dihedral_style before AngleAngleTorsion Coeffs");
"Must define dihedral_style before "
"AngleAngleTorsion Coeffs");
if (firstpass) dihedralcoeffs(4); if (firstpass) dihedralcoeffs(4);
else skip_lines(ndihedraltypes); else skip_lines(ndihedraltypes);
} else if (strcmp(keyword,"BondBond13 Coeffs") == 0) { } else if (strcmp(keyword,"BondBond13 Coeffs") == 0) {
if (atom->avec->dihedrals_allow == 0) if (atom->avec->dihedrals_allow == 0)
error->all(FLERR,"Invalid data file section: BondBond13 Coeffs"); error->all(FLERR,"Invalid data file section: BondBond13 Coeffs");
if (force->dihedral == nullptr) if (force->dihedral == nullptr)
error->all(FLERR, error->all(FLERR,"Must define dihedral_style before BondBond13 Coeffs");
"Must define dihedral_style before BondBond13 Coeffs");
if (firstpass) dihedralcoeffs(5); if (firstpass) dihedralcoeffs(5);
else skip_lines(ndihedraltypes); else skip_lines(ndihedraltypes);
@ -712,8 +728,7 @@ void ReadData::command(int narg, char **arg)
if (atom->avec->impropers_allow == 0) if (atom->avec->impropers_allow == 0)
error->all(FLERR,"Invalid data file section: AngleAngle Coeffs"); error->all(FLERR,"Invalid data file section: AngleAngle Coeffs");
if (force->improper == nullptr) if (force->improper == nullptr)
error->all(FLERR, error->all(FLERR,"Must define improper_style before AngleAngle Coeffs");
"Must define improper_style before AngleAngle Coeffs");
if (firstpass) impropercoeffs(1); if (firstpass) impropercoeffs(1);
else skip_lines(nimpropertypes); else skip_lines(nimpropertypes);
@ -724,8 +739,7 @@ void ReadData::command(int narg, char **arg)
for (i = 0; i < nfix; i++) for (i = 0; i < nfix; i++)
if (strcmp(keyword,fix_section[i]) == 0) { if (strcmp(keyword,fix_section[i]) == 0) {
if (firstpass) fix(fix_index[i],keyword); if (firstpass) fix(fix_index[i],keyword);
else skip_lines(modify->fix[fix_index[i]]-> else skip_lines(modify->fix[fix_index[i]]->read_data_skip_lines(keyword));
read_data_skip_lines(keyword));
parse_keyword(0); parse_keyword(0);
break; break;
} }
@ -882,8 +896,7 @@ void ReadData::command(int narg, char **arg)
bigint nblocal = atom->nlocal; bigint nblocal = atom->nlocal;
MPI_Allreduce(&nblocal,&natoms,1,MPI_LMP_BIGINT,MPI_SUM,world); MPI_Allreduce(&nblocal,&natoms,1,MPI_LMP_BIGINT,MPI_SUM,world);
if (natoms != atom->natoms) if (natoms != atom->natoms)
error->all(FLERR, error->all(FLERR,"Read_data shrink wrap did not assign all atoms correctly");
"Read_data shrink wrap did not assign all atoms correctly");
} }
// restore old styles, when reading with nocoeff flag given // restore old styles, when reading with nocoeff flag given
@ -947,17 +960,6 @@ void ReadData::header(int firstpass)
atom->nimpropertypes = extra_improper_types; atom->nimpropertypes = extra_improper_types;
} }
// customize for new sections
const char *section_keywords[NSECTIONS] =
{"Atoms","Velocities","Ellipsoids","Lines","Triangles","Bodies",
"Bonds","Angles","Dihedrals","Impropers",
"Masses","Pair Coeffs","PairIJ Coeffs","Bond Coeffs","Angle Coeffs",
"Dihedral Coeffs","Improper Coeffs",
"BondBond Coeffs","BondAngle Coeffs","MiddleBondTorsion Coeffs",
"EndBondTorsion Coeffs","AngleTorsion Coeffs",
"AngleAngleTorsion Coeffs","BondBond13 Coeffs","AngleAngle Coeffs"};
// skip 1st line of file // skip 1st line of file
if (me == 0) { if (me == 0) {
@ -1186,9 +1188,7 @@ void ReadData::header(int firstpass)
// check that exiting string is a valid section keyword // check that exiting string is a valid section keyword
parse_keyword(1); parse_keyword(1);
for (n = 0; n < NSECTIONS; n++) if (!is_data_section(keyword))
if (strcmp(keyword,section_keywords[n]) == 0) break;
if (n == NSECTIONS)
error->all(FLERR,"Unknown identifier in data file: {}",keyword); error->all(FLERR,"Unknown identifier in data file: {}",keyword);
// error checks on header values // error checks on header values
@ -1730,8 +1730,7 @@ void ReadData::bodies(int firstpass, AtomVec *ptr)
error->one(FLERR,"Too many values in body lines in data file"); error->one(FLERR,"Too many values in body lines in data file");
if (onebody+1 > MAXBODY) if (onebody+1 > MAXBODY)
error->one(FLERR, error->one(FLERR,"Too many lines in one body in data file - boost MAXBODY");
"Too many lines in one body in data file - boost MAXBODY");
nchunk++; nchunk++;
nline += onebody+1; nline += onebody+1;

View File

@ -29,6 +29,7 @@ class ReadData : public Command {
ReadData(class LAMMPS *); ReadData(class LAMMPS *);
~ReadData(); ~ReadData();
void command(int, char **); void command(int, char **);
static bool is_data_section(const std::string &);
private: private:
int me, compressed; int me, compressed;