[Rdkit-discuss] editing molecules

2014-03-01 Thread S.L. Chan
Good evening,

I have some molecules with neutral amidine heads. I would like to
make them positively charged. So I do:

mol = Chem.MolFromMolFile('input.mol', removeHs=False)
q = Chem.MolFromSmarts('[NH]=C[NH2]')
for mat in mol.GetSubstructMatches(q):
    mol.GetAtomWithIdx(mat[0]).SetFormalCharge(1)
    mol.GetAtomWithIdx(mat[0]).SetNumExplicitHs(2)

m1 = Chem.AddHs(mol, addCoords=True)

However, RDKit complains that the valence for the N in concern is 4,
which is greater than permitted.

I also tried it without the "SetNumExplicitHs" line. While it did not
crash, it refused to add any H to the charged N.

Attached is an example input file. But I got the same result with
other input files too.

Have I missed some steps?

Thank you.

Ling


1LPZ.mol
Description: Binary data
--
Flow-based real-time traffic analytics software. Cisco certified tool.
Monitor traffic, SLAs, QoS, Medianet, WAAS etc. with NetFlow Analyzer
Customize your own dashboards, set traffic alerts and generate reports.
Network behavioral analysis & security monitoring. All-in-one tool.
http://pubads.g.doubleclick.net/gampad/clk?id=126839071&iu=/4140/ostg.clktrk___
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss


Re: [Rdkit-discuss] editing molecules

2014-03-03 Thread Greg Landrum
Hi,


On Sun, Mar 2, 2014 at 3:57 AM, S.L. Chan  wrote:

> Good evening,
>
> I have some molecules with neutral amidine heads. I would like to
> make them positively charged. So I do:
>
> mol = Chem.MolFromMolFile('input.mol', removeHs=False)
> q = Chem.MolFromSmarts('[NH]=C[NH2]')
> for mat in mol.GetSubstructMatches(q):
> mol.GetAtomWithIdx(mat[0]).SetFormalCharge(1)
> mol.GetAtomWithIdx(mat[0]).SetNumExplicitHs(2)
>
> m1 = Chem.AddHs(mol, addCoords=True)
>
> However, RDKit complains that the valence for the N in concern is 4,
> which is greater than permitted.
>

You need to sanitize the molecule after you change the charge. This will
automatically adjust the H count for you:

In [5]: mol = Chem.MolFromSmiles('CC(=N)N')
In [6]: q = Chem.MolFromSmarts('[NH]=C[NH2]')
In [7]: for mat in mol.GetSubstructMatches(q):
mol.GetAtomWithIdx(mat[0]).SetFormalCharge(1)
   ...:

In [8]: Chem.SanitizeMol(mol)
Out[8]: rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_NONE

In [11]: mh=Chem.AddHs(mol)

In [12]: print Chem.MolToSmiles(mol)
CC(N)=[NH2+]

In [13]: print Chem.MolToSmiles(mh)
[H]N([H])C(=[N+]([H])[H])C([H])([H])[H]


Best,
-greg





>
> I also tried it without the "SetNumExplicitHs" line. While it did not
> crash, it refused to add any H to the charged N.
>
> Attached is an example input file. But I got the same result with
> other input files too.
>
> Have I missed some steps?
>
> Thank you.
>
> Ling
>
>
>
> --
> Flow-based real-time traffic analytics software. Cisco certified tool.
> Monitor traffic, SLAs, QoS, Medianet, WAAS etc. with NetFlow Analyzer
> Customize your own dashboards, set traffic alerts and generate reports.
> Network behavioral analysis & security monitoring. All-in-one tool.
>
> http://pubads.g.doubleclick.net/gampad/clk?id=126839071&iu=/4140/ostg.clktrk
> ___
> Rdkit-discuss mailing list
> Rdkit-discuss@lists.sourceforge.net
> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
>
--
Subversion Kills Productivity. Get off Subversion & Make the Move to Perforce.
With Perforce, you get hassle-free workflows. Merge that actually works. 
Faster operations. Version large binaries.  Built-in WAN optimization and the
freedom to use Git, Perforce or both. Make the move to Perforce.
http://pubads.g.doubleclick.net/gampad/clk?id=122218951&iu=/4140/ostg.clktrk___
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss