Thanks, that's worked brilliantly. Seems quicker than my original approach too.

On 21/01/2016 15:17, Brian Kelley wrote:
Without a concrete example, this solution may not be appropriate, but I believe the function you want is "ReplaceCore".

ReplaceCore(...)

ReplaceCore( (Mol)mol, (Mol)coreQuery [, (bool)replaceDummies=True [, (bool)labelByIndex=False [, (bool)requireDummyMatch=False]]]) -> Mol :

Removes the core of a molecule and labels the sidechains with dummy atoms.



I just have python available currently so this may not be appropriate, but here goes:

>>> m1 = Chem.MolFromSmiles("Cc1ccccc1N")

>>> m2 = Chem.MolFromSmiles("c1ccccc1")

>>> mcs = MCS.FindMCS([m1, m2])

>>> frag = Chem.ReplaceCore(m1, Chem.MolFromSmarts(mcs.smarts))

>>> print "SideChains:", Chem.MolToSmiles(frag)

SideChains: [*]C.[*]N


I hope this helps (at least the steps).

Now if you are just trying to extract side chains from the results of reactions, we have recently added helper functions to solve that (They should be exposed in the next release).


ReduceProductToSideChains(...)

ReduceProductToSideChains( (Mol)product [, (bool)addDummyAtoms=True]) -> Mol :

reduce the product of a reaction to the side chains added by the reaction.

The output is a molecule with attached wildcards indicating where the product was attached. The isotope of the dummy atom is the reaction map number of the product's atom (if available).


If this would be useful, let us know, I would be happy to have a tester prior to release.

Brian Kelley

On Thu, Jan 21, 2016 at 9:41 AM, James Wallace <chp11...@sheffield.ac.uk <mailto:chp11...@sheffield.ac.uk>> wrote:

    Hi,
    I'm using the KNIME implementation to write my own nodes, and I'm
    running into an issue. For the process I'm trying to do I'm trying to
    subtract the MCS between two molecules away from the larger
    molecule, to
    leave a list of fragments. I'm aware of the substructure matching, but
    I'm not sure how to subtract the matching atoms from a molecule graph
    within RDKit. As I say, I'm working with the Java version, but any
    pointers towards the fucntions needed would be useful. At the moment
    I've got (in pseudo code)

             RWMol mol1a = RWMol.MolFromSmiles(reactant_string, 0, true);
             RWMol mol2a = RWMol.MolFromSmiles(product_string, 0, true);

         frag_bonds = mol2a.GetSubstructMatches(mol1a);

    But I'm unsure as to what to do with the array of matches to achieve
    what I want. Can I strip out the dummy atoms automatically, or is this
    something that is best achieved by processing the SMILES string?

    
------------------------------------------------------------------------------
    Site24x7 APM Insight: Get Deep Visibility into Application Performance
    APM + Mobile APM + RUM: Monitor 3 App instances at just $35/Month
    Monitor end-to-end web transactions and take corrective actions now
    Troubleshoot faster and improve end-user experience. Signup Now!
    http://pubads.g.doubleclick.net/gampad/clk?id=267308311&iu=/4140
    _______________________________________________
    Rdkit-discuss mailing list
    Rdkit-discuss@lists.sourceforge.net
    <mailto:Rdkit-discuss@lists.sourceforge.net>
    https://lists.sourceforge.net/lists/listinfo/rdkit-discuss



------------------------------------------------------------------------------
Site24x7 APM Insight: Get Deep Visibility into Application Performance
APM + Mobile APM + RUM: Monitor 3 App instances at just $35/Month
Monitor end-to-end web transactions and take corrective actions now
Troubleshoot faster and improve end-user experience. Signup Now!
http://pubads.g.doubleclick.net/gampad/clk?id=267308311&iu=/4140
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to