BioSimSpace.Align.matchAtoms#

BioSimSpace.Align.matchAtoms(molecule0, molecule1, scoring_function='rmsd_align', matches=1, return_scores=False, prematch={}, timeout=5.0000 secs, complete_rings_only=True, max_scoring_matches=1000, property_map0={}, property_map1={})[source]#

Find mappings between atom indices in molecule0 to those in molecule1. Molecules are aligned using a Maximum Common Substructure (MCS) search. When requesting more than one match, the mappings will be sorted using a scoring function and returned in order of best to worst score. (Note that, depending on the scoring function the “best” score may have the lowest value.).

Parameters:
  • molecule0 (Molecule) – The molecule of interest.

  • molecule1 (Molecule) – The reference molecule.

  • scoring_function (str) –

    The scoring function used to match atoms. Available options are:
    • ”rmsd”

      Calculate the root mean squared distance between the coordinates of atoms in molecule0 to those that they map to in molecule1.

    • ”rmsd_align”

      Align molecule0 to molecule1 based on the mapping before computing the above RMSD score.

    • ”rmsd_flex_align”

      Flexibly align molecule0 to molecule1 based on the mapping before computing the above RMSD score. (based on the ‘fkcombu’. package: https://pdbj.org/kcombu)

  • matches (int) – The maximum number of matches to return. (Sorted in order of score).

  • return_scores (bool) – Whether to return a list containing the scores for each mapping.

  • prematch (dict) – A dictionary of atom mappings that must be included in the match.

  • timeout (BioSimSpace.Types.Time) – The timeout for the maximum common substructure search.

  • complete_rings_only (bool) – Whether to only match complete rings during the MCS search. This option is only relevant to MCS performed using RDKit and will be ignored when falling back on Sire.

  • max_scoring_matches (int) – The maximum number of matching MCS substructures to consider when computing mapping scores. Consider reducing this if you find the matchAtoms function to be taking an excessive amount of time. This option is only relevant to MCS performed using RDKit and will be ignored when falling back on Sire.

  • property_map0 (dict) – A dictionary that maps “properties” in molecule0 to their user defined values. This allows the user to refer to properties with their own naming scheme, e.g. { “charge” : “my-charge” }

  • property_map1 (dict) – A dictionary that maps “properties” in molecule1 to their user defined values.

Returns:

matches – The best atom mapping, a list containing a user specified number of the best mappings ranked by their score, or a tuple containing the list of best mappings and a list of the corresponding scores.

Return type:

dict, [dict], ([dict], list)

Examples

Find the best maximum common substructure mapping between two molecules.

>>> import BioSimSpace as BSS
>>> mapping = BSS.Align.matchAtoms(molecule0, molecule1)

Find the 5 best mappings.

>>> import BioSimSpace as BSS
>>> mappings = BSS.Align.matchAtoms(molecule0, molecule1, matches=5)

Find the 5 best mappings along with their ranking scores.

>>> import BioSimSpace as BSS
>>> mappings, scores = BSS.Align.matchAtoms(molecule0, molecule1, matches=5, return_scores=True)

Find the 5 best mappings along with their ranking scores. Score by flexibly aligning molecule0 to molecule1 based on each mapping and computing the root mean squared displacement of the matched atoms.

>>> import BioSimSpace as BSS
>>> mappings, scores = BSS.Align.matchAtoms(molecule0, molecule1, matches=5, return_scores=True, scoring_function="rmsd_flex_align")

Find the best mapping that contains a prematch (this is a dictionary mapping atom indices in molecule0 to those in molecule1).

>>> import BioSimSpace as BSS
>>> mapping = BSS.Align.matchAtoms(molecule0, molecule1, prematch={0 : 10, 3 : 7})