On Oct 31, 2019, at 12:49, LJ <lijo89j...@gmail.com> wrote:
> How to find out the terminal ring (ring having no connections with any other
> ring) using SMARTS notation?

That is not possible with SMARTS.

For one, SMARTS matches a fixed substructure (possibly recursively), meaning 
there is no way to match a ring of indefinite size. You can only match a ring 
of a given, fixed size.

For another, you cannot express "no connections with" in SMARTS, only 
"connection with something which isn't"

For example, the SMARTS "[R1]@1@[R1]@[R1]@[R1]@[R1]@[R1]1~[!R1]" will match a 
6-membered ring where the atoms are only in a single SSSR ring (not in a fused 
ring), and where one of the atoms is connected with something which isn't in a 
ring.

However, consider "C1CCCCC1(O)c1ccccc1" where the cyclohexane is connected to 
an hydroxyl and to a benzo.

The above SMARTS matches, even though there are two rings connected to each 
other.

Here's one way to approach what you are looking for.

Open Babel has a way to find the smallest set of smallest rings ("SSSR"). All 
disconnected rings must be in the SSSR.

The "[R1]" SMARTS pattern will find all atoms which are in only a single SSSR 
ring.

The pattern "[R]!@[R]" will find all ring atoms which are connected to another 
ring atom but not through a ring bond.

So, you want the SSSR rings which contain only atoms which are in a single SSSR 
ring, and which are not bonded to another ring atom through a non-ring bond.

The following does something like what I think you want:


import openbabel
if openbabel.__version__.startswith("2."):
    import pybel
else:
    from openbabel import pybel # Use the Open Babel 3.0 API

ring_atoms_pat = pybel.Smarts("[R1]")
connected_ring_atoms_pat = pybel.Smarts("[R]!@[R]")

def find_disconnected_rings(mol):
    rings = mol.sssr

    # Get the atoms which are only in a single SSSR ring
    single_ring_atoms = set(m[0] for m in ring_atoms_pat.findall(mol))
    
    # Find all ring atoms which are connected to another ring
    # atom through a non-ring bond.
    connected_ring_atoms = set()
    for m in connected_ring_atoms_pat.findall(mol):
        connected_ring_atoms.add(m[0])
        connected_ring_atoms.add(m[1])
 
    # Get the rings which contain atoms which are only in a single
    # SSSR ring and which are not connected to another ring.
    disconnected_rings = []
    for ring in rings:
        if all(atom in single_ring_atoms and atom not in connected_ring_atoms
                     for atom in ring._path):
            disconnected_rings.append(ring)
    return disconnected_rings

for smiles in (
        "c1ccccc1O",
        "CNP",
        "c1ccccc1c1ccccc1",
        "C1CCCCC1(O)c1ccccc1",
        "C1CCCCC1(O)Pc1ccccc1",
        "C1CCCCC1(O)Pc1ccccc1-C1CCC1",
        ):
    mol = pybel.readstring("smi", smiles)
    print("SMILES:", smiles)
    has_ring = False
    for ring_atoms in find_disconnected_rings(mol):
        print("  ", ring_atoms._path)
        has_ring = True
    if not has_ring:
        print("  (no disconnected rings)")

The output of this is:


SMILES: c1ccccc1O
   (3, 4, 5, 6, 1, 2)
SMILES: CNP
  (no disconnected rings)
SMILES: c1ccccc1c1ccccc1
  (no disconnected rings)
SMILES: C1CCCCC1(O)c1ccccc1
  (no disconnected rings)
SMILES: C1CCCCC1(O)Pc1ccccc1
   (3, 4, 5, 6, 1, 2)
   (11, 12, 13, 14, 9, 10)
SMILES: C1CCCCC1(O)Pc1ccccc1-C1CCC1
   (3, 4, 5, 6, 1, 2)


                                Andrew
                                da...@dalkescientific.com




_______________________________________________
OpenBabel-discuss mailing list
OpenBabel-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/openbabel-discuss

Reply via email to