Hi Paolo,
Thanks, this is an interesting feature - may come useful one day!
I found a different solution - via torsion constraints:

mp = ChemicalForceFields.MMFFGetMoleculeProperties(mol)
ff = ChemicalForceFields.MMFFGetMoleculeForceField(mol, mp)
sp=Chem.MolFromSmarts("c1ccsc1!@c") #thiophene connected to another
aromatic ring

atoms=(1,0,4,5) #out-of-plane bend, atoms 1, 0 and 4 are in the thiophene
ring, atom 5 is the first atom of the next ring

maplist = mol.GetSubstructMatches(sp)

if (len(maplist)>0):
        for match in maplist :

            angle=180.0 #desired dihedral angle
            a=[]
            for i in range (0,4) :
                a.append(match[atoms[i]])
            cur_angle=rdMolTransforms.GetDihedralDeg(mol.GetConformer(),
a[0],a[1],a[2], a[3])

            if abs(abs(cur_angle)-abs(angle))>45.0: # closer to 0 than to
180?
                angle=0.0 # freeze at 0
   ff.MMFFAddTorsionConstraint(a[0],a[1],a[2],a[3], False, angle, angle,
1e4)

ff.Minimize()


Best wishes,

Michal

On 15 April 2015 at 23:06, Paolo Tosco <paolo.to...@unito.it> wrote:

>  Dear Michal,
>
> please find attached a small script which accomplishes what you describe
> by a different approach, i.e. it minimizes only the methyl group in
> 2-methylthiophene while keeping the rest fixed, effectively pushing it back
> in plane. Would that work for you?
>
> Best,
> Paolo
>
>
> On 04/15/2015 10:57 AM, Michal Krompiec wrote:
>
>  Hello,
> I'm trying to manipulate out-of-plane bends of substituents attached to
> aromatic rings. Obviously, they should lie in the plane of the ring, but
> sometimes (after constrained optimization) they come slightly out of plane
> and there doesn't seem to be a way to push/rotate them back to the ring
> plane. Adding constraints to the forcefield (which seems to work for
> out-of-plane bends as well) is not a perfect solution.
>
>  For example, there is no way to manipulate the out-of-plane bend
> (dihedral angle) of the methyl group in the 2-methylthiophene molecule
> c1cc(C)sc1, defined by the first four atoms c1cc(C). If you try to use
> rdMolTransforms.SetDihedralDeg, it will raise an exception because the bond
> around it tries to rotate is part of a ring (line 375 in MolTransforms.cpp)
> - and for a good reason, because what I wanted to do is to rotate just the
> methyl group (atom l as defined in setDihedralRad, and everything attached
> to it), not all atoms connected to atom k being on the "right hand side" of
> the j-k bond.
>
>  Would it make sense to make a modified version of setDihedralRad():
> setOutOfPlaneBendRad(), with the following modifications:
> * omit line 375: if(queryIsBondInRing(bondJK)) throw ValueErrorException("bond
> (j,k) must not belong to a ring"); instead perhaps check if bond k,l
> belongs to a ring
> * change line 401: _toBeMovedIdxList(mol, jAtomId, kAtomId, alist); to
> _toBeMovedIdxList(mol, kAtomId, lAtomId, alist);  (we want to rotate the
> tree rooted in atom l, not atom k)
>
>  Best regards,
> Michal
>
>
> ------------------------------------------------------------------------------
> BPM Camp - Free Virtual Workshop May 6th at 10am PDT/1PM EDT
> Develop your own process in accordance with the BPMN 2 standard
> Learn Process modeling best practices with Bonita BPM through live 
> exerciseshttp://www.bonitasoft.com/be-part-of-it/events/bpm-camp-virtual- 
> event?utm_
> source=Sourceforge_BPM_Camp_5_6_15&utm_medium=email&utm_campaign=VA_SF
>
>
>
> _______________________________________________
> Rdkit-discuss mailing 
> listRdkit-discuss@lists.sourceforge.nethttps://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
>
>
>
> ------------------------------------------------------------------------------
> BPM Camp - Free Virtual Workshop May 6th at 10am PDT/1PM EDT
> Develop your own process in accordance with the BPMN 2 standard
> Learn Process modeling best practices with Bonita BPM through live
> exercises
> http://www.bonitasoft.com/be-part-of-it/events/bpm-camp-virtual-
> event?utm_
> source=Sourceforge_BPM_Camp_5_6_15&utm_medium=email&utm_campaign=VA_SF
> _______________________________________________
> Rdkit-discuss mailing list
> Rdkit-discuss@lists.sourceforge.net
> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
>
------------------------------------------------------------------------------
BPM Camp - Free Virtual Workshop May 6th at 10am PDT/1PM EDT
Develop your own process in accordance with the BPMN 2 standard
Learn Process modeling best practices with Bonita BPM through live exercises
http://www.bonitasoft.com/be-part-of-it/events/bpm-camp-virtual- event?utm_
source=Sourceforge_BPM_Camp_5_6_15&utm_medium=email&utm_campaign=VA_SF
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to