Skip to content

Commit

Permalink
Fixes #5923 (#7039)
Browse files Browse the repository at this point in the history
  • Loading branch information
greglandrum committed Jan 19, 2024
1 parent 8127c85 commit 7ba3be0
Show file tree
Hide file tree
Showing 4 changed files with 108 additions and 9 deletions.
46 changes: 46 additions & 0 deletions Code/GraphMol/MolInterchange/molinterchange_catch.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -377,3 +377,49 @@ TEST_CASE("Test segv reported in #6890") {
FileParseException);
}
}

TEST_CASE("github #5923: add more error checking to substance groups") {
SECTION("basics") {
std::string moljson =
R"JSON({"rdkitjson":{"version":11},"defaults":{"atom":{"z":6,"impHs":0,"chg":0,"nRad":0,"isotope":0,"stereo":"unspecified"},"bond":{"bo":1,"stereo":"unspecified"}},
"molecules":[{"name":"","atoms":[{"z":0},{"impHs":1},{"z":8},{"z":0},{"impHs":3}],
"bonds":[{"atoms":[0,1]},{"atoms":[1,2]},{"atoms":[2,3]},{"atoms":[1,4]}],
"substanceGroups":[{"properties":{"TYPE":"SRU","index":1,"CONNECT":"HT","LABEL":"n","DATAFIELDS":"[]"},
"atoms":[2,1,4],"bonds":[2,0],"parentAtoms":[0,3],
"brackets":[[[-3.9538,4.3256,0.0],[-3.0298,2.7252,0.0],[0.0,0.0,0.0]],[[-5.4618,2.8611,0.0],[-6.3858,4.4615,0.0],[0.0,0.0,0.0]]]}],
"conformers":[{"dim":2,"coords":[[-6.7083,3.2083],[-5.3747,3.9783],[-4.041,3.2083],[-2.7073,3.9783],[-5.3747,5.5183]]}],"extensions":[{"name":"rdkitRepresentation","formatVersion":2,"toolkitVersion":"2023.09.4","cipRanks":[0,3,4,1,2]},{"name":"rdkitQueries","formatVersion":10,"toolkitVersion":"2023.09.4","atomQueries":[{"descr":"AtomNull","tag":40},{},{},{"descr":"AtomNull","tag":40},{}]}]}]})JSON";
// everything ok
{
auto mols = MolInterchange::JSONDataToMols(moljson);
REQUIRE(mols.size() == 1);
CHECK(getSubstanceGroups(*mols[0]).size() == 1);
}
// bad atom index
{
std::string badjson = moljson;
std::string lookFor = R"JSON("atoms":[2,1,4])JSON";
badjson.replace(badjson.find(lookFor), lookFor.size(),
R"JSON("atoms":[2,1,9])JSON");
CHECK_THROWS_AS(MolInterchange::JSONDataToMols(badjson),
ValueErrorException);
}
// bad parentAtom index
{
std::string badjson = moljson;
std::string lookFor = R"JSON("parentAtoms":[0,3])JSON";
badjson.replace(badjson.find(lookFor), lookFor.size(),
R"JSON("parentAtoms":[8,3])JSON");
CHECK_THROWS_AS(MolInterchange::JSONDataToMols(badjson),
ValueErrorException);
}
// bad bond index
{
std::string badjson = moljson;
std::string lookFor = R"JSON("bonds":[2,0])JSON";
badjson.replace(badjson.find(lookFor), lookFor.size(),
R"JSON("bonds":[2,7])JSON");
CHECK_THROWS_AS(MolInterchange::JSONDataToMols(badjson),
ValueErrorException);
}
}
}
50 changes: 46 additions & 4 deletions Code/GraphMol/SubstanceGroup.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -60,10 +60,48 @@ unsigned int SubstanceGroup::getIndexInMol() const {
return sgroupItr - sgroups.begin();
}

void SubstanceGroup::addAtomWithIdx(unsigned int idx) {
void SubstanceGroup::setAtoms(std::vector<unsigned int> atoms) {
PRECONDITION(dp_mol, "bad mol");
auto natoms = dp_mol->getNumAtoms();
for (auto atomIdx : atoms) {
if (atomIdx >= natoms) {
std::ostringstream errout;
errout << "Atom index " << atomIdx << " is out of range";
throw ValueErrorException(errout.str());
}
}
d_atoms = std::move(atoms);
}
void SubstanceGroup::setParentAtoms(std::vector<unsigned int> patoms) {
PRECONDITION(dp_mol, "bad mol");
auto natoms = dp_mol->getNumAtoms();
for (auto atomIdx : patoms) {
if (atomIdx >= natoms) {
std::ostringstream errout;
errout << "Atom index " << atomIdx << " is out of range";
throw ValueErrorException(errout.str());
}
}
d_patoms = std::move(patoms);
}
void SubstanceGroup::setBonds(std::vector<unsigned int> bonds) {
PRECONDITION(dp_mol, "bad mol");
PRECONDITION(dp_mol->getAtomWithIdx(idx), "wrong atom index");
auto nbonds = dp_mol->getNumBonds();
for (auto bondIdx : bonds) {
if (bondIdx >= nbonds) {
std::ostringstream errout;
errout << "Bond index " << bondIdx << " is out of range";
throw ValueErrorException(errout.str());
}
}
d_bonds = std::move(bonds);
}

void SubstanceGroup::addAtomWithIdx(unsigned int idx) {
PRECONDITION(dp_mol, "bad mol");
if (idx >= dp_mol->getNumAtoms()) {
throw ValueErrorException("Atom index out of range");
}
d_atoms.push_back(idx);
}

Expand All @@ -76,7 +114,9 @@ void SubstanceGroup::addAtomWithBookmark(int mark) {

void SubstanceGroup::addParentAtomWithIdx(unsigned int idx) {
PRECONDITION(dp_mol, "bad mol");

if (idx >= dp_mol->getNumAtoms()) {
throw ValueErrorException("Atom index out of range");
}
if (std::find(d_atoms.begin(), d_atoms.end(), idx) == d_atoms.end()) {
std::ostringstream errout;
errout << "Atom " << idx << " is not a member of current SubstanceGroup";
Expand All @@ -103,7 +143,9 @@ void SubstanceGroup::addParentAtomWithBookmark(int mark) {

void SubstanceGroup::addBondWithIdx(unsigned int idx) {
PRECONDITION(dp_mol, "bad mol");
PRECONDITION(dp_mol->getBondWithIdx(idx), "wrong bond index");
if (idx >= dp_mol->getNumBonds()) {
throw ValueErrorException("Bond index out of range");
}

d_bonds.push_back(idx);
}
Expand Down
8 changes: 3 additions & 5 deletions Code/GraphMol/SubstanceGroup.h
Original file line number Diff line number Diff line change
Expand Up @@ -171,11 +171,9 @@ class RDKIT_GRAPHMOL_EXPORT SubstanceGroup : public RDProps {
const std::vector<unsigned int> &getParentAtoms() const { return d_patoms; }
const std::vector<unsigned int> &getBonds() const { return d_bonds; }

void setAtoms(std::vector<unsigned int> atoms) { d_atoms = std::move(atoms); }
void setParentAtoms(std::vector<unsigned int> patoms) {
d_patoms = std::move(patoms);
}
void setBonds(std::vector<unsigned int> bonds) { d_bonds = std::move(bonds); }
void setAtoms(std::vector<unsigned int> atoms);
void setParentAtoms(std::vector<unsigned int> patoms);
void setBonds(std::vector<unsigned int> bonds);

const std::vector<Bracket> &getBrackets() const { return d_brackets; }
const std::vector<CState> &getCStates() const { return d_cstates; }
Expand Down
13 changes: 13 additions & 0 deletions Code/GraphMol/catch_sgroups.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,7 @@
#include <GraphMol/Chirality.h>
#include <GraphMol/SmilesParse/SmilesParse.h>
#include "GraphMol/FileParsers/FileParsers.h"
#include "GraphMol/MolInterchange/MolInterchange.h"

using namespace RDKit;

Expand Down Expand Up @@ -507,4 +508,16 @@ M END
TEST_ASSERT(index == count);
}
}
}

TEST_CASE("github #5923: add more error checking to substance groups") {
auto mol = "CCCO"_smiles;
REQUIRE(mol);
SubstanceGroup sg(mol.get(), "SUP");
CHECK_THROWS_AS(sg.addAtomWithIdx(4), ValueErrorException);
CHECK_THROWS_AS(sg.addParentAtomWithIdx(4), ValueErrorException);
CHECK_THROWS_AS(sg.addBondWithIdx(4), ValueErrorException);
CHECK_THROWS_AS(sg.setAtoms({1, 4}), ValueErrorException);
CHECK_THROWS_AS(sg.setParentAtoms({1, 4}), ValueErrorException);
CHECK_THROWS_AS(sg.setBonds({1, 4}), ValueErrorException);
}

0 comments on commit 7ba3be0

Please sign in to comment.