bond/react: update error messages

This commit is contained in:
jrgissing 2019-05-30 23:27:23 -06:00
parent e3e5a962b0
commit 80d906d445
2 changed files with 72 additions and 55 deletions

View File

@ -308,8 +308,7 @@ FixBondReact::FixBondReact(LAMMPS *lmp, int narg, char **arg) :
onemol->check_attributes(0);
twomol->check_attributes(0);
if (onemol->natoms != twomol->natoms)
error->all(FLERR,"Post-reacted template must contain the same "
"number of atoms as the pre-reacted template");
error->all(FLERR,"Bond/react: Reaction templates must contain the same number of atoms");
get_molxspecials();
read(i);
fclose(fp);
@ -324,7 +323,7 @@ FixBondReact::FixBondReact(LAMMPS *lmp, int narg, char **arg) :
delete [] files;
if (atom->molecular != 1)
error->all(FLERR,"Cannot use fix bond/react with non-molecular systems");
error->all(FLERR,"Bond/react: Cannot use fix bond/react with non-molecular systems");
// check if bonding atoms are 1-2, 1-3, or 1-4 bonded neighbors
// if so, we don't need non-bonded neighbor list
@ -665,7 +664,7 @@ void FixBondReact::init()
// check cutoff for iatomtype,jatomtype
for (int i = 0; i < nreacts; i++) {
if (force->pair == NULL || cutsq[i][1] > force->pair->cutsq[iatomtype[i]][jatomtype[i]])
error->all(FLERR,"Fix bond/react cutoff is longer than pairwise cutoff");
error->all(FLERR,"Bond/react: Fix bond/react cutoff is longer than pairwise cutoff");
}
// need a half neighbor list, built every Nevery steps
@ -1174,7 +1173,7 @@ void FixBondReact::superimpose_algorithm()
// let's go ahead and catch the simplest of hangs
//if (hang_catch > onemol->natoms*4)
if (hang_catch > atom->nlocal*30) {
error->one(FLERR,"Excessive iteration of superimpose algorithm");
error->one(FLERR,"Bond/react: Excessive iteration of superimpose algorithm");
}
}
}
@ -1287,7 +1286,7 @@ void FixBondReact::make_a_guess()
for (int i = 0; i < nxspecial[atom->map(glove[pion][1])][0]; i++) {
if (atom->map(xspecial[atom->map(glove[pion][1])][i]) < 0) {
error->all(FLERR,"Fix bond/react needs ghost atoms from further away1"); // parallel issues.
error->all(FLERR,"Bond/react: Fix bond/react needs ghost atoms from further away1"); // parallel issues.
}
if (i_limit_tags[(int)atom->map(xspecial[atom->map(glove[pion][1])][i])] != 0) {
status = GUESSFAIL;
@ -1398,7 +1397,7 @@ void FixBondReact::check_a_neighbor()
//another check for ghost atoms. perhaps remove the one in make_a_guess
if (atom->map(glove[(int)onemol_xspecial[pion][neigh]-1][1]) < 0) {
error->all(FLERR,"Fix bond/react needs ghost atoms from further away2");
error->all(FLERR,"Bond/react: Fix bond/react needs ghost atoms from further away2");
}
for (int j = 0; j < onemol_nxspecial[onemol_xspecial[pion][neigh]-1][0]; j++) {
@ -1450,7 +1449,7 @@ void FixBondReact::check_a_neighbor()
//another check for ghost atoms. perhaps remove the one in make_a_guess
if (atom->map(glove[(int)onemol_xspecial[pion][neigh]-1][1]) < 0) {
error->all(FLERR,"Fix bond/react needs ghost atoms from further away3");
error->all(FLERR,"Bond/react: Fix bond/react needs ghost atoms from further away3");
}
for (int ii = 0; ii < onemol_nxspecial[onemol_xspecial[pion][neigh]-1][0]; ii++) {
@ -1492,7 +1491,7 @@ void FixBondReact::crosscheck_the_neighbor()
glove[onemol_xspecial[pion][trace]-1][0] == 0) {
if (avail_guesses == MAXGUESS) {
error->warning(FLERR,"Fix bond/react failed because MAXGUESS set too small. ask developer for info");
error->warning(FLERR,"Bond/react: Fix bond/react failed because MAXGUESS set too small. ask developer for info");
status = GUESSFAIL;
return;
}
@ -1561,7 +1560,7 @@ void FixBondReact::inner_crosscheck_loop()
//another check for ghost atoms. perhaps remove the one in make_a_guess
if (atom->map(glove[(int)onemol_xspecial[pion][neigh]-1][1]) < 0) {
error->all(FLERR,"Fix bond/react needs ghost atoms from further away4");
error->all(FLERR,"Bond/react: Fix bond/react needs ghost atoms from further away4");
}
if (guess_branch[avail_guesses-1] == 0) avail_guesses--;
@ -1722,7 +1721,7 @@ void FixBondReact::find_landlocked_atoms(int myrxn)
// if atoms change types, but aren't landlocked, that's bad
for (int i = 0; i < twomol->natoms; i++) {
if (twomol->type[i] != onemol->type[equivalences[i][1][myrxn]-1] && landlocked_atoms[i][myrxn] == 0)
error->one(FLERR,"Atom affected by reaction too close to template edge");
error->one(FLERR,"Bond/react: Atom affected by reaction too close to template edge");
}
// additionally, if a bond changes type, but neither involved atom is landlocked, bad
@ -1738,7 +1737,7 @@ void FixBondReact::find_landlocked_atoms(int myrxn)
onemol_batom = onemol->bond_atom[onemol_atomi-1][m];
if (onemol_batom == equivalences[twomol_atomj-1][1][myrxn]) {
if (twomol->bond_type[i][j] != onemol->bond_type[onemol_atomi-1][m]) {
error->one(FLERR,"Bond type affected by reaction too close to template edge");
error->one(FLERR,"Bond/react: Bond type affected by reaction too close to template edge");
}
}
}
@ -1748,7 +1747,7 @@ void FixBondReact::find_landlocked_atoms(int myrxn)
onemol_batom = onemol->bond_atom[onemol_atomj-1][m];
if (onemol_batom == equivalences[i][1][myrxn]) {
if (twomol->bond_type[i][j] != onemol->bond_type[onemol_atomj-1][m]) {
error->one(FLERR,"Bond type affected by reaction too close to template edge");
error->one(FLERR,"Bond/react: Bond type affected by reaction too close to template edge");
}
}
}
@ -1764,7 +1763,7 @@ void FixBondReact::find_landlocked_atoms(int myrxn)
int ii = reverse_equiv[i][1][myrxn] - 1;
for (int j = 0; j < twomol_nxspecial[ii][0]; j++) {
if (delete_atoms[equivalences[twomol_xspecial[ii][j]-1][1][myrxn]-1][myrxn] == 0) {
error->one(FLERR,"A deleted atom cannot remain bonded to an atom that is not deleted");
error->one(FLERR,"Bond/react: A deleted atom cannot remain bonded to an atom that is not deleted");
}
}
}
@ -1775,7 +1774,7 @@ void FixBondReact::find_landlocked_atoms(int myrxn)
for (int i = 0; i < twomol->natoms; i++) {
if (twomol_nxspecial[i][0] != onemol_nxspecial[equivalences[i][1][myrxn]-1][0] && landlocked_atoms[i][myrxn] == 0) {
char str[128];
sprintf(str,"An atom in 'react #%d' changes bond connectivity but not atom type",myrxn+1);
sprintf(str,"Bond/react: An atom in 'react #%d' changes bond connectivity but not atom type",myrxn+1);
error->warning(FLERR,str);
break;
}
@ -2263,7 +2262,7 @@ void FixBondReact::update_everything()
if (landlocked_atoms[j][rxnID] == 1) {
for (int k = 0; k < nspecial[atom->map(update_mega_glove[jj+1][i])][2]; k++) {
if (atom->map(special[atom->map(update_mega_glove[jj+1][i])][k]) < 0) {
error->all(FLERR,"Fix bond/react needs ghost atoms from further away - most likely too many processors");
error->all(FLERR,"Bond/react: Fix bond/react needs ghost atoms from further away - most likely too many processors");
}
}
}
@ -2388,7 +2387,7 @@ void FixBondReact::update_everything()
bond_atom[atom->map(update_mega_glove[jj+1][i])][insert_num] = update_mega_glove[equivalences[twomol->bond_atom[j][p]-1][1][rxnID]][i];
num_bond[atom->map(update_mega_glove[jj+1][i])]++;
if (num_bond[atom->map(update_mega_glove[jj+1][i])] > atom->bond_per_atom)
error->one(FLERR,"Bond/react bonds/atom exceed system bonds/atom");
error->one(FLERR,"Bond/react topology/atom exceed system topology/atom");
delta_bonds++;
}
}
@ -2464,7 +2463,7 @@ void FixBondReact::update_everything()
angle_atom3[atom->map(update_mega_glove[jj+1][i])][insert_num] = update_mega_glove[equivalences[twomol->angle_atom3[j][p]-1][1][rxnID]][i];
num_angle[atom->map(update_mega_glove[jj+1][i])]++;
if (num_angle[atom->map(update_mega_glove[jj+1][i])] > atom->angle_per_atom)
error->one(FLERR,"Bond/react angles/atom exceed system angles/atom");
error->one(FLERR,"Bond/react topology/atom exceed system topology/atom");
delta_angle++;
}
}
@ -2547,7 +2546,7 @@ void FixBondReact::update_everything()
dihedral_atom4[atom->map(update_mega_glove[jj+1][i])][insert_num] = update_mega_glove[equivalences[twomol->dihedral_atom4[j][p]-1][1][rxnID]][i];
num_dihedral[atom->map(update_mega_glove[jj+1][i])]++;
if (num_dihedral[atom->map(update_mega_glove[jj+1][i])] > atom->dihedral_per_atom)
error->one(FLERR,"Bond/react dihedrals/atom exceed system dihedrals/atom");
error->one(FLERR,"Bond/react topology/atom exceed system topology/atom");
delta_dihed++;
}
}
@ -2630,7 +2629,7 @@ void FixBondReact::update_everything()
improper_atom4[atom->map(update_mega_glove[jj+1][i])][insert_num] = update_mega_glove[equivalences[twomol->improper_atom4[j][p]-1][1][rxnID]][i];
num_improper[atom->map(update_mega_glove[jj+1][i])]++;
if (num_improper[atom->map(update_mega_glove[jj+1][i])] > atom->improper_per_atom)
error->one(FLERR,"Bond/react impropers/atom exceed system impropers/atom");
error->one(FLERR,"Bond/react topology/atom exceed system topology/atom");
delta_imprp++;
}
}
@ -2728,7 +2727,7 @@ void FixBondReact::read(int myrxn)
// skip 1st line of file
eof = fgets(line,MAXLINE,fp);
if (eof == NULL) error->one(FLERR,"Unexpected end of superimpose file");
if (eof == NULL) error->one(FLERR,"Bond/react: Unexpected end of superimpose file");
// read header lines
// skip blank lines or lines that start with "#"
@ -2782,7 +2781,7 @@ void FixBondReact::read(int myrxn)
DeleteAtoms(line, myrxn);
} else if (strcmp(keyword,"Constraints") == 0) {
Constraints(line, myrxn);
} else error->one(FLERR,"Unknown section in superimpose file");
} else error->one(FLERR,"Bond/react: Unknown section in map file");
parse_keyword(1,line,keyword);
@ -2790,13 +2789,13 @@ void FixBondReact::read(int myrxn)
// error check
if (bondflag == 0 || equivflag == 0)
error->all(FLERR,"Superimpose file missing BondingIDs or Equivalences section\n");
error->all(FLERR,"Bond/react: Map file missing BondingIDs or Equivalences section\n");
if (update_edges_flag[myrxn] == 2 && customedgesflag == 0)
error->all(FLERR,"Map file must have a Custom Edges section when using 'update_edges custom'\n");
error->all(FLERR,"Bond/react: Map file must have a Custom Edges section when using 'update_edges custom'\n");
if (update_edges_flag[myrxn] != 2 && customedgesflag == 1)
error->all(FLERR,"Specify 'update_edges custom' to include Custom Edges section in map file\n");
error->all(FLERR,"Bond/react: Specify 'update_edges custom' to include Custom Edges section in map file\n");
}
void FixBondReact::EdgeIDs(char *line, int myrxn)
@ -2842,7 +2841,7 @@ void FixBondReact::CustomEdges(char *line, int myrxn)
else if (strcmp(edgemode,"charges") == 0)
custom_edges[tmp-1][myrxn] = 1;
else
error->one(FLERR,"Illegal value in 'Custom Edges' section of map file");
error->one(FLERR,"Bond/react: Illegal value in 'Custom Edges' section of map file");
}
delete [] edgemode;
}
@ -2873,7 +2872,7 @@ void FixBondReact::Constraints(char *line, int myrxn)
constraints[myrxn][3] = tmp[2]*tmp[2]; // using square of distance
constraints[myrxn][4] = tmp[3]*tmp[3];
} else
error->one(FLERR,"Illegal constraint type in 'Constraints' section of map file");
error->one(FLERR,"Bond/react: Illegal constraint type in 'Constraints' section of map file");
}
delete [] constraint_type;
}
@ -2883,7 +2882,7 @@ void FixBondReact::open(char *file)
fp = fopen(file,"r");
if (fp == NULL) {
char str[128];
snprintf(str,128,"Cannot open superimpose file %s",file);
snprintf(str,128,"Bond/react: Cannot open map file %s",file);
error->one(FLERR,str);
}
}
@ -2896,7 +2895,7 @@ void FixBondReact::readline(char *line)
else n = strlen(line) + 1;
}
MPI_Bcast(&n,1,MPI_INT,0,world);
if (n == 0) error->all(FLERR,"Unexpected end of superimpose file");
if (n == 0) error->all(FLERR,"Bond/react: Unexpected end of map file");
MPI_Bcast(line,n,MPI_CHAR,0,world);
}

View File

@ -189,47 +189,65 @@ Self-explanatory. Check the input script syntax and compare to the
documentation for the command. You can use -echo screen as a
command-line option when running LAMMPS to see the offending line.
E: Invalid exclude group name
E: Bond/react: Cannot use fix bond/react with non-molecular systems
Exclude group name should not previously be defined.
E: Cannot use fix bond/react with non-molecular systems
Only systems with bonds that can be changed can be used. Atom_style
Only systems with bonds that can be changed can be used. Atom_style
template does not qualify.
E: Fix bond/react cutoff is longer than pairwise cutoff
E: Bond/react: Rmax cutoff is longer than pairwise cutoff
This is not allowed because bond creation is done using the
pairwise neighbor list.
This is not allowed because bond creation is done using the pairwise
neighbor list.
E: Molecule template ID for fix bond/react does not exist
E: Bond/react: Molecule template ID for fix bond/react does not exist
A valid molecule template must have been created with the molecule command.
A valid molecule template must have been created with the molecule
command.
E: Superimpose file errors:
E: Bond/react: Reaction templates must contain the same number of atoms
Please ensure superimpose file is properly formatted.
There should be a one-to-one correspondence between atoms in the
pre-reacted and post-reacted templates, as specified by the map file.
E: Atom affected by reaction too close to template edge
E: Bond/react: Unknown section in map file
Please ensure reaction map files are properly formatted.
E: Bond/react: Atom affected by reaction too close to template edge
This means an atom which changes type during the reaction is too close
to an 'edge' atom defined in the superimpose file. This could cause incorrect
assignment of bonds, angle, etc. Generally, this means you must include
more atoms in your templates, such that there are at least two atoms
between each atom involved in the reaction and an edge atom.
to an 'edge' atom defined in the superimpose file. This could cause
incorrect assignment of bonds, angle, etc. Generally, this means you
must include more atoms in your templates, such that there are at
least two atoms between each atom involved in the reaction and an edge
atom.
E: Fix bond/react needs ghost atoms from farther away
E: Bond/react: Fix bond/react needs ghost atoms from farther away
This is because a processor needs to superimpose the entire unreacted
molecule template onto simulation atoms it can 'see.' The comm_modify cutoff
command can be used to extend the communication range.
molecule template onto simulation atoms it knows about. The
comm_modify cutoff command can be used to extend the communication
range.
E: Excessive iteration of superimpose algorithm
E: Bond/react: A deleted atom cannot remain bonded to an atom that is not deleted
You may have discovered a bug! But first, please double check that your
molecule template atom types, bond types, etc. are consistent with your simulation,
and that all atoms affected by a reaction are sufficently separated from edge atoms.
If this issue persists, please contact the developer.
Self-explanatory.
W: Bond/react: An atom in 'react #%d' changes bond connectivity but not atom type
You may want to double-check that all atom types are properly assigned
in the post-reaction template.
E: Bond/react special bond generation overflow
The number of special bonds per-atom created by a reaction exceeds the
system setting. See the read_data or create_box command for how to
specify this value.
E: Bond/react topology/atom exceed system topology/atom
The number of bonds, angles etc per-atom created by a reaction exceeds
the system setting. See the read_data or create_box command for how to
specify this value.
*/