Dear Paolo Thank you for this illustrative example. I tried the same method on a few other molecules and didn't find consistent results. In the notebook attached, the first molecule is left with a bivalent carbon, and the second molecule is only sanitized after explicitly calling the SanitizeFrags option when getting individual fragments. Am I missing something?
Kind regards Puck On Tue, 18 Feb 2020 at 11:17, Paolo Tosco <paolo.tosco.m...@gmail.com> wrote: > Hi Puck, > > I modified your previous example to use FragmentOnBonds(): > > from rdkit import Chemfrom rdkit.Chem.Draw import IPythonConsole, > MolsToGridImage > > smiles = "COC"mol = Chem.MolFromSmiles("COC")mol > > Here there are no hydrogens (as expected from the SMILES) > > mol = Chem.AddHs(mol)mol > > mol_f = Chem.FragmentOnBonds(mol, (2,)) > > mol_f > > for a in mol_f.GetAtoms(): > if (a.GetAtomicNum() == 0): > # convert the dummy into hydrogen > a.SetAtomicNum(1) > # if you do not need the labels, you may clear them > a.SetIsotope(0) > > mol_f > > Cheers, > p. > > On 18/02/2020 09:45, Puck van Gerwen wrote: > > Dear Paolo > > Thanks very much for your continuous help. I'm not sure I like the idea of > explicitly defining the valence, so I'll consider switching to > FragmentOnBonds(). From the example you gave previously though, it quenched > the fragments with dummy atoms rather than hydrogen. Is it possible to > quench with hydrogen using FragmentOnBonds()? > > Kind regards > Puck > > On Mon, 17 Feb 2020 at 13:42, Paolo Tosco <paolo.tosco.m...@gmail.com> > wrote: > >> Hi Puck, >> >> sorry dor the delay in replying. >> >> After removing the bond you will need to adjust the number of explicit Hs >> on both ends with SetNumExplicitHs(), and then add those Hs in the >> molecule graph as real atoms with Chem.AddHs(): >> >> rwmol = Chem.RWMol(mol) >> bonds = rwmol.GetBonds() >> b = rwmol.GetBondWithIdx(2) >> a1 = b.GetBeginAtom() >> a2 = b.GetEndAtom() >> rwmol.RemoveBond(a1.GetIdx(), a2.GetIdx()) >> a1.SetNumExplicitHs(a1.GetNumExplicitHs() + 1) >> a2.SetNumExplicitHs(a2.GetNumExplicitHs() + 1) >> rwmol = Chem.AddHs(rwmol, explicitOnly=True) >> >> The rdmolops.FragmentOnBonds() functions adjusts the valences for you, >> but you will need to adjust them manually if you do the fragmentation at a >> lower level with RemoveBond(). >> >> I have attached the modified notebook. >> >> Cheers, >> p. >> On 12/02/2020 10:03, Puck van Gerwen wrote: >> >> Dear all, >> >> I am trying to read in SMILES to generate mol objects which I then break >> into fragments using rwmol.RemoveBond(). Thereafter I want to sanitize the >> fragments by saturating with hydrogen. However, I am finding that rdkit >> often doesn't sanitize the fragments consistently, leaving trivalent carbon >> atoms. I've attached a jupyter notebook of an example. Could anyone help me >> consistently generate sanitized fragments? >> >> Best regards >> -- >> *Puck van Gerwen* >> Doktorandin / PhD candidate >> Departement Chemie >> Klingelbergstrasse 80 >> CH-4056 Basel >> https://chemspacelab.org/ >> >> >> _______________________________________________ >> Rdkit-discuss mailing >> listRdkit-discuss@lists.sourceforge.nethttps://lists.sourceforge.net/lists/listinfo/rdkit-discuss >> >> > > -- > *Puck van Gerwen* > Doktorandin / PhD candidate > Departement Chemie > Klingelbergstrasse 80 > CH-4056 Basel > https://chemspacelab.org/ > > -- *Puck van Gerwen* Doktorandin / PhD candidate Departement Chemie Klingelbergstrasse 80 CH-4056 Basel https://chemspacelab.org/
TestFragmentation.ipynb
Description: Binary data
_______________________________________________ Rdkit-discuss mailing list Rdkit-discuss@lists.sourceforge.net https://lists.sourceforge.net/lists/listinfo/rdkit-discuss