RDKit
Open-source cheminformatics and machine learning.
Loading...
Searching...
No Matches
EnumerateStereoisomers.h
Go to the documentation of this file.
1//
2// Copyright (C) David Cosgrove 2025.
3//
4// @@ All Rights Reserved @@
5// This file is part of the RDKit.
6// The contents are covered by the terms of the BSD license
7// which is included in the file license.txt, found at the root
8// of the RDKit source tree.
9//
10// This is and the corresponding C++ file contain a transliteration
11// of EnumerateStereoisomers.py. It enumerates all tetrahedral,
12// double bond and astropisomer stereochemistry. The last being
13// and enhancement over the Python original.
14
15#ifndef RD_ENUMERATESTEREOISOMERS_H
16#define RD_ENUMERATESTEREOISOMERS_H
17
19
20#include <random>
21#include <unordered_set>
22
23#include <boost/dynamic_bitset/dynamic_bitset.hpp>
24
25#include <RDGeneral/export.h>
26#include <GraphMol/ROMol.h>
27#include <GraphMol/RWMol.h>
28
29namespace RDKit {
31
33 bool tryEmbedding{false}; // If true, the process attempts to generate
34 // a standard RDKit distance geometry
35 // conformation for the stereoisomer. If this
36 // fails, we assume that the stereoisomer is
37 // non-physical and don't return it. NOTE that
38 // this is computationally expensive and is just
39 // a heuristic that could result in
40 // stereoisomers being lost.
41 bool onlyUnassigned{true}; // If true, stereocenters which have a
42 // specified stereochemistry will not be
43 // perturbed unless they are part of a relative
44 // stereo group.
45 bool onlyStereoGroups{false}; // If true, only find stereoisomers that
46 // differ at the StereoGroups associated
47 // with the molecule.
48 bool unique{true}; // If true, only stereoisomers that differ in canonical
49 // SMILES will be returned.
50 std::uint64_t maxIsomers{0}; // The maximum number of isomers to yield.
51 // If the number of possible isomers is
52 // greater than maxIsomers, a random subset
53 // will be yielded. If 0, all isomers are
54 // yielded. Since every additional
55 // stereocenter doubles the number of results
56 // (and execution time) it's important to
57 // keep an eye on this.
58 int randomSeed{-1}; // Seed for random number generator. -1 means don't
59 // seed.
60};
61
62// Class that enumerates the stereoisomers of a molecule. Acts like a
63// Python generator so in principle has no limit on the number of stereoisomers
64// it can produce.
66 public:
69 const ROMol &mol,
71 bool verbose = false);
76 delete;
78
79 std::uint64_t getStereoisomerCount() const;
80
81 // Return another stereoisomer, or an empty unique_ptr if we're done.
82 std::unique_ptr<ROMol> next();
83
84 private:
85 RWMol d_mol;
86 const StereoEnumerationOptions d_options;
87 bool d_verbose;
88 // The number of isomers successfully returned so far. This may be
89 // lower than d_seen.size() if d_options.tryEmbedding is true and there
90 // have been embedding failures.
91 std::uint64_t d_numReturned{0};
92 // The number we need to return, the smaller of 2**N, where N is the
93 // number of stereo centers, or d_options.maxIsomers.
94 std::uint64_t d_numToReturn{1024};
95 // 2**N.
96 std::uint64_t d_totalPoss{0};
97 std::unordered_set<std::string> d_generatedIsomers;
98
99 // For the random bools
100 std::unique_ptr<std::mt19937> d_randGen;
101 std::bernoulli_distribution d_randDis{0.5};
102
103 // Classes for setting the orientation at a particular stereocenter.
104 std::vector<std::unique_ptr<details::Flipper>> d_flippers;
105
106 // The stereo orientations we've already made
107 std::unordered_set<boost::dynamic_bitset<>> d_seen;
108
109 void buildFlippers();
110 std::unique_ptr<ROMol> generateRandomIsomer();
111 bool embeddable(ROMol &isomer);
112};
113} // namespace EnumerateStereoisomers
114} // namespace RDKit
115
116#endif // RD_ENUMERATESTEREOISOMERS_H
Defines the primary molecule class ROMol as well as associated typedefs.
Defines the editable molecule class RWMol.
StereoisomerEnumerator & operator=(StereoisomerEnumerator &&other)=delete
StereoisomerEnumerator(const ROMol &mol, const StereoEnumerationOptions &options=StereoEnumerationOptions(), bool verbose=false)
StereoisomerEnumerator(StereoisomerEnumerator &&other)=delete
StereoisomerEnumerator(const StereoisomerEnumerator &other)=delete
StereoisomerEnumerator & operator=(const StereoisomerEnumerator &other)=delete
RWMol is a molecule class that is intended to be edited.
Definition RWMol.h:32
#define RDKIT_ENUMERATESTEREOISOMERS_EXPORT
Definition export.h:161
Std stuff.