On Fri, Feb 8, 2019 at 5:52 AM Noki Lee <noki.le...@gmail.com> wrote:

> Thanks a lot!
>
> It's working for my project!
>

glad to hear it.


> But, in the later way, STEREOTRANS (STEREOCIS) can be replaced by STEREOE
> (STEREOZ), can't it?
> I mean it might cause a problem that we consider one of N atom becomes C
> atom or both.
>

You should always be fine using STEREOTRANS and STEREOCIS. The problem with
STEREOE and STEREOZ is that you (technically) need to know the relative CIP
priorities of the stereo atoms; this problem is why we added STEREOTRANS
and STEREOCIS.

-greg



> Thanks, greg
>
> On Thu, Feb 7, 2019 at 2:31 PM Greg Landrum <greg.land...@gmail.com>
> wrote:
>
>> Hi,
>>
>> This has gotten a little bit easier in the most recent RDKit version.
>> Here's a simple example:
>>
>> In [3]: print(rdkit.__version__)
>> 2018.09.1
>> In [4]: from rdkit import Chem
>> In [5]: m = Chem.MolFromSmiles('CN=NC')
>> In [6]: m.GetBondWithIdx(1).SetStereoAtoms(0,3)
>> In [7]: m.GetBondWithIdx(1).SetStereo(Chem.BondStereo.STEREOCIS)
>> In [8]: Chem.AssignStereochemistry(m,force=True)
>> In [9]: Chem.MolToSmiles(m)
>> Out[9]: 'C/N=N\\C'
>>
>>
>> Or a more complicated one that's closer to what it looks like you want to
>> do:
>>
>> In [10]: m = Chem.MolFromSmiles('c1ccccc1N=Nc1ccccc1')
>> In [11]: match = m.GetSubstructMatch(Chem.MolFromSmarts('cN=Nc'))
>> In [12]: b = m.GetBondBetweenAtoms(match[1],match[2])
>> In [13]: b.SetStereoAtoms(match[0],match[3])
>> In [14]: b.SetStereo(Chem.BondStereo.STEREOTRANS)
>> In [15]: Chem.AssignStereochemistry(m,force=True)
>> In [16]: Chem.MolToSmiles(m)
>> Out[16]: 'c1ccc(/N=N/c2ccccc2)cc1'
>>
>>
>> I hope that helps,
>> -greg
>>
>>
>>
>>
>>
>> On Wed, Feb 6, 2019 at 12:31 PM Noki Lee <noki.le...@gmail.com> wrote:
>>
>>> Hi, RDkit.
>>>
>>> I'm looking for a way to change cis-trans for a molecule at the diimide
>>> part.
>>> There was a case for the Alkene E-Z change.
>>>
>>> https://sourceforge.net/p/rdkit/mailman/message/35011276/
>>>
>>> But when I apply above method, it doesn't work for the diimide case.
>>>
>>> For examples, let say Mymol be made from azobenzene:
>>> c1ccc(N=Nc2ccccc2)cc1
>>> I wrote down the procedure:
>>> 1. Detect substructure of N=N part by GetSubstructMatch
>>> 2. And then apply ENDDOWNRIGHT and ENDUPRIGHT for each atomIdx
>>> but the result is same for the original one.
>>>
>>> Thanks all!
>>>
>>> _______________________________________________
>>> Rdkit-discuss mailing list
>>> Rdkit-discuss@lists.sourceforge.net
>>> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>>>
>>
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to