Hi Nick (and Paolo), It's not the molecules that trigger the infinite loop (I tried it with one of mine and it hanged as well). The problem is that no coordinates are added for the hydrogens in the script. If you set addCoords=True when adding the hydrogens, the minimization works (at least it did with me...).
for mol in suppl:
molList.append(Chem.AddHs(mol, addCoords=True))
Best,
Sereina
2013/10/8 Paolo Tosco <[email protected]>
> Hi Nick,
>
> would you mind sending me the SD file which is triggering the infinite
> loop? Then I'll come back to you as soon as I find out something.
>
> Best,
> p.
>
>
>
> On 10/08/2013 11:36 AM, Nicholas Firth wrote:
>
> Hi RDKitters,
>
> I'm having an issue using the MMFF to minimise a CORINA conformation.
> I've written a little script which adds hydrogens to a molecule then tries
> to use the MMFF forcefield to minimise the conformer. The problem is that
> the script hangs on the minimise step.
>
> This error only occurs when I add hydrogens to the conformation, I
> assume the reason for this is because the hydrogens are all added at the
> origin. Is there a way of getting round this (in RDKit, as I want to keep
> the AddHs function)?
>
> I've included the script below, it will work like this however if you
> switch wrong the commented lines in the first for loop then the it no
> longer works.
>
>
> from rdkit import Chem
> from rdkit.Chem import AllChem
> from rdkit.Chem import ChemicalForceFields
> from sys import argv
>
> suppl = Chem.SDMolSupplier(argv[1])
> molList = []
>
> for mol in suppl:
> #molList.append(Chem.AddHs(mol))
> molList.append(mol)
> del suppl
>
> #w = Chem.SDWriter(argv[1])
> w = Chem.SDWriter('test.sdf')
>
> for mol in molList:
> mp = ChemicalForceFields.MMFFGetMoleculeProperties(mol)
> field = AllChem.MMFFGetMoleculeForceField(mol, mp)
> field.Minimize()
> w.write(mol)
>
> w.close()
>
>
> Thanks in advance.
>
> Best,
> Nick
>
> *Nicholas C. Firth* | PhD Student | Cancer Therapeutics****
> The Institute of Cancer Research | 15 Cotswold Road | Belmont | Sutton |
> Surrey | SM2 5NG****
>
> *T* 020 8722 4033 | *E* [email protected] | *W* www.icr.ac.uk | *
> Twitter* @ICRnews <https://twitter.com/ICRnews>****
>
> *Facebook* www.facebook.com/theinstituteofcancerresearch****
>
> *Making the discoveries that defeat cancer*
>
>
>
> The Institute of Cancer Research: Royal Cancer Hospital, a charitable
> Company Limited by Guarantee, Registered in England under Company No.
> 534147 with its Registered Office at 123 Old Brompton Road, London SW7 3RP.
>
> This e-mail message is confidential and for use by the addressee only. If
> the message is received by anyone other than the addressee, please return
> the message to the sender by replying to it and then delete the message
> from your computer and network.
>
>
> ------------------------------------------------------------------------------
> October Webinars: Code for Performance
> Free Intel webinars can help you accelerate application performance.
> Explore tips for MPI, OpenMP, advanced profiling, and more. Get the most from
> the latest Intel processors and coprocessors. See abstracts and register
> >http://pubads.g.doubleclick.net/gampad/clk?id=60134071&iu=/4140/ostg.clktrk
>
>
>
> _______________________________________________
> Rdkit-discuss mailing
> [email protected]https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
>
>
> --
> ==========================================================
> Paolo Tosco, Ph.D.
> Department of Drug Science and Technology
> Via Pietro Giuria, 9 - 10125 Torino (Italy)
> Tel: +39 011 670 7680 | Mob: +39 348 5537206
> Fax: +39 011 670 7687 | E-mail: [email protected]http://open3dqsar.org |
> http://open3dalign.org
> ==========================================================
>
>
>
> ------------------------------------------------------------------------------
> October Webinars: Code for Performance
> Free Intel webinars can help you accelerate application performance.
> Explore tips for MPI, OpenMP, advanced profiling, and more. Get the most
> from
> the latest Intel processors and coprocessors. See abstracts and register >
> http://pubads.g.doubleclick.net/gampad/clk?id=60134071&iu=/4140/ostg.clktrk
> _______________________________________________
> Rdkit-discuss mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
>
<<image/gif>>
------------------------------------------------------------------------------ October Webinars: Code for Performance Free Intel webinars can help you accelerate application performance. Explore tips for MPI, OpenMP, advanced profiling, and more. Get the most from the latest Intel processors and coprocessors. See abstracts and register > http://pubads.g.doubleclick.net/gampad/clk?id=60134071&iu=/4140/ostg.clktrk
_______________________________________________ Rdkit-discuss mailing list [email protected] https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

