File: MHFPWrapper.cpp

package info (click to toggle)
rdkit 202503.1-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 220,160 kB
  • sloc: cpp: 399,240; python: 77,453; ansic: 25,517; java: 8,173; javascript: 4,005; sql: 2,389; yacc: 1,565; lex: 1,263; cs: 1,081; makefile: 580; xml: 229; fortran: 183; sh: 105
file content (240 lines) | stat: -rw-r--r-- 11,577 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
//
//  2019, Daniel Probst, Reymond Group @ University of Bern
//
//   @@ All Rights Reserved @@
//  This file is part of the RDKit.
//  The contents are covered by the terms of the BSD license
//  which is included in the file license.txt, found at the root
//  of the RDKit source tree.
//

#include <boost/python.hpp>
#include <RDBoost/Wrap.h>
#include <GraphMol/Fingerprints/MHFP.h>
#include <vector>

namespace python = boost::python;
using RDKit::MHFPFingerprints::MHFPEncoder;

namespace RDKit {
namespace MHFPWrapper {

typedef std::vector<std::vector<uint32_t>> VectMinHashVect;
typedef std::vector<ExplicitBitVect> VectExplicitBitVect;

template <typename T>
std::vector<T> ListToVector(const python::object &obj) {
  return std::vector<T>(python::stl_input_iterator<T>(obj),
                        python::stl_input_iterator<T>());
}

UINT_VECT
FromStringArray(MHFPEncoder *mhfpEnc, python::list &vec) {
  std::vector<std::string> vec_tmp = ListToVector<std::string>(vec);
  return mhfpEnc->FromStringArray(vec_tmp);
}

UINT_VECT
FromArray(MHFPEncoder *mhfpEnc, python::list &vec) {
  std::vector<uint32_t> vec_tmp = ListToVector<uint32_t>(vec);
  return mhfpEnc->FromArray(vec_tmp);
}

STR_VECT
CreateShinglingFromSmiles(MHFPEncoder *mhfpEnc, std::string smiles,
                          unsigned char radius = 3, bool rings = true,
                          bool isomeric = false, bool kekulize = true,
                          unsigned char min_radius = 1) {
  return mhfpEnc->CreateShingling(smiles, radius, rings, isomeric, kekulize,
                                  min_radius);
}

STR_VECT
CreateShinglingFromMol(MHFPEncoder *mhfpEnc, ROMol mol,
                       unsigned char radius = 3, bool rings = true,
                       bool isomeric = false, bool kekulize = true,
                       unsigned char min_radius = 1) {
  return mhfpEnc->CreateShingling(mol, radius, rings, isomeric, kekulize,
                                  min_radius);
}

UINT_VECT
EncodeSmiles(MHFPEncoder *mhfpEnc, std::string smiles, unsigned char radius = 3,
             bool rings = true, bool isomeric = false, bool kekulize = true,
             unsigned char min_radius = 1) {
  return mhfpEnc->Encode(smiles, radius, rings, isomeric, kekulize, min_radius);
}

UINT_VECT
EncodeMol(MHFPEncoder *mhfpEnc, ROMol mol, unsigned char radius = 3,
          bool rings = true, bool isomeric = false, bool kekulize = true,
          unsigned char min_radius = 1) {
  return mhfpEnc->Encode(mol, radius, rings, isomeric, kekulize, min_radius);
}

VectMinHashVect EncodeSmilesBulk(MHFPEncoder *mhfpEnc, python::list &smiles,
                                 unsigned char radius = 3, bool rings = true,
                                 bool isomeric = false, bool kekulize = true,
                                 unsigned char min_radius = 1) {
  std::vector<std::string> vec = ListToVector<std::string>(smiles);
  return mhfpEnc->Encode(vec, radius, rings, isomeric, kekulize, min_radius);
}

// There are access problems for vector_indexing_suite and std::vector<ROMol>.
// So let's fallback to a unefficient Python list.
VectMinHashVect EncodeMolsBulk(MHFPEncoder *mhfpEnc, python::list &mols,
                               unsigned char radius = 3, bool rings = true,
                               bool isomeric = false, bool kekulize = true,
                               unsigned char min_radius = 1) {
  std::vector<ROMol> vec = ListToVector<ROMol>(mols);
  return mhfpEnc->Encode(vec, radius, rings, isomeric, kekulize, min_radius);
}

// SECFP

ExplicitBitVect EncodeSECFPSmiles(MHFPEncoder *mhfpEnc, std::string smiles,
                                  unsigned char radius = 3, bool rings = true,
                                  bool isomeric = false, bool kekulize = true,
                                  unsigned char min_radius = 1,
                                  size_t length = 2048) {
  return mhfpEnc->EncodeSECFP(smiles, radius, rings, isomeric, kekulize,
                              min_radius, length);
}

ExplicitBitVect EncodeSECFPMol(MHFPEncoder *mhfpEnc, ROMol mol,
                               unsigned char radius = 3, bool rings = true,
                               bool isomeric = false, bool kekulize = true,
                               unsigned char min_radius = 1,
                               size_t length = 2048) {
  return mhfpEnc->EncodeSECFP(mol, radius, rings, isomeric, kekulize,
                              min_radius, length);
}

VectExplicitBitVect EncodeSECFPSmilesBulk(
    MHFPEncoder *mhfpEnc, python::list &smiles, unsigned char radius = 3,
    bool rings = true, bool isomeric = false, bool kekulize = true,
    unsigned char min_radius = 1, size_t length = 2048) {
  std::vector<std::string> vec = ListToVector<std::string>(smiles);
  return mhfpEnc->EncodeSECFP(vec, radius, rings, isomeric, kekulize,
                              min_radius, length);
}

VectExplicitBitVect EncodeSECFPMolsBulk(
    MHFPEncoder *mhfpEnc, python::list &mols, unsigned char radius = 3,
    bool rings = true, bool isomeric = false, bool kekulize = true,
    unsigned char min_radius = 1, size_t length = 2048) {
  std::vector<ROMol> vec = ListToVector<ROMol>(mols);
  return mhfpEnc->EncodeSECFP(vec, radius, rings, isomeric, kekulize,
                              min_radius, length);
}

BOOST_PYTHON_FUNCTION_OVERLOADS(CreateShinglingFromSmilesOverloads,
                                CreateShinglingFromSmiles, 2, 7)
BOOST_PYTHON_FUNCTION_OVERLOADS(CreateShinglingFromMolOverloads,
                                CreateShinglingFromMol, 2, 7)

BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeSmilesOverloads, EncodeSmiles, 2, 7)
BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeMolOverloads, EncodeMol, 2, 7)
BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeSmilesBulkOverloads, EncodeSmilesBulk, 2,
                                7)
BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeMolsBulkOverloads, EncodeMolsBulk, 2, 7)

BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeSECFPSmilesOverloads, EncodeSECFPSmiles,
                                2, 8)
BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeSECFPMolOverloads, EncodeSECFPMol, 2, 8)
BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeSECFPSmilesBulkOverloads,
                                EncodeSECFPSmilesBulk, 2, 8)
BOOST_PYTHON_FUNCTION_OVERLOADS(EncodeSECFPMolsBulkOverloads,
                                EncodeSECFPMolsBulk, 2, 8)

BOOST_PYTHON_MODULE(rdMHFPFingerprint) {
  python::class_<MHFPEncoder>(
      "MHFPEncoder", python::init<python::optional<unsigned int, unsigned int>>(
                         python::args("self", "n_permutations", "seed")))
      .def("FromStringArray", FromStringArray, python::args("self", "vec"),
           "Creates a MHFP vector from a list of arbitrary strings.")
      .def("FromArray", FromArray, python::args("self", "vec"),
           "Creates a MHFP vector from a list of unsigned integers.")
      .def("CreateShinglingFromSmiles", CreateShinglingFromSmiles,
           CreateShinglingFromSmilesOverloads(
               (python::arg("self"), python::arg("smiles"),
                python::arg("radius") = 3, python::arg("rings") = true,
                python::arg("isomeric") = false,
                python::arg("kekulize") = false, python::arg("min_radius") = 1),
               "Creates a shingling (a list of circular n-grams / "
               "substructures) from a SMILES string."))
      .def("CreateShinglingFromMol", CreateShinglingFromMol,
           CreateShinglingFromMolOverloads(
               (python::arg("self"), python::arg("mol"),
                python::arg("radius") = 3, python::arg("rings") = true,
                python::arg("isomeric") = false,
                python::arg("kekulize") = false, python::arg("min_radius") = 1),
               "Creates a shingling (a list of circular n-grams / "
               "substructures) from a RDKit Mol instance."))
      .def("EncodeSmiles", EncodeSmiles,
           EncodeSmilesOverloads(
               (python::arg("self"), python::arg("smiles"),
                python::arg("radius") = 3, python::arg("rings") = true,
                python::arg("isomeric") = false,
                python::arg("kekulize") = false, python::arg("min_radius") = 1),
               "Creates a MHFP vector from a SMILES string."))
      .def("EncodeMol", EncodeMol,
           EncodeMolOverloads(
               (python::arg("self"), python::arg("mol"),
                python::arg("radius") = 3, python::arg("rings") = true,
                python::arg("isomeric") = false,
                python::arg("kekulize") = false, python::arg("min_radius") = 1),
               "Creates a MHFP vector from an RDKit Mol instance."))
      .def("EncodeSmilesBulk", EncodeSmilesBulk,
           EncodeSmilesBulkOverloads(
               (python::arg("self"), python::arg("smiles"),
                python::arg("radius") = 3, python::arg("rings") = true,
                python::arg("isomeric") = false,
                python::arg("kekulize") = false, python::arg("min_radius") = 1),
               "Creates a MHFP vector from a list of SMILES strings."))
      .def("EncodeMolsBulk", EncodeMolsBulk,
           EncodeMolsBulkOverloads(
               (python::arg("self"), python::arg("mols"),
                python::arg("radius") = 3, python::arg("rings") = true,
                python::arg("isomeric") = false,
                python::arg("kekulize") = false, python::arg("min_radius") = 1),
               "Creates a MHFP vector from a list of RDKit Mol instances."))
      .def(
          "EncodeSECFPSmiles", EncodeSECFPSmiles,
          EncodeSECFPSmilesOverloads(
              (python::arg("self"), python::arg("smiles"),
               python::arg("radius") = 3, python::arg("rings") = true,
               python::arg("isomeric") = false, python::arg("kekulize") = false,
               python::arg("min_radius") = 1, python::arg("length") = 2048),
              "Creates a SECFP binary vector from a SMILES string."))
      .def(
          "EncodeSECFPMol", EncodeSECFPMol,
          EncodeSECFPMolOverloads(
              (python::arg("self"), python::arg("smiles"),
               python::arg("radius") = 3, python::arg("rings") = true,
               python::arg("isomeric") = false, python::arg("kekulize") = false,
               python::arg("min_radius") = 1, python::arg("length") = 2048),
              "Creates a SECFP binary vector from an RDKit Mol instance."))
      .def(
          "EncodeSECFPSmilesBulk", EncodeSECFPSmilesBulk,
          EncodeSECFPSmilesBulkOverloads(
              (python::arg("self"), python::arg("smiles"),
               python::arg("radius") = 3, python::arg("rings") = true,
               python::arg("isomeric") = false, python::arg("kekulize") = false,
               python::arg("min_radius") = 1, python::arg("length") = 2048),
              "Creates a SECFP binary vector from a list of SMILES strings."))
      .def(
          "EncodeSECFPMolsBulk", EncodeSECFPMolsBulk,
          EncodeSECFPMolsBulkOverloads(
              (python::arg("self"), python::arg("smiles"),
               python::arg("radius") = 3, python::arg("rings") = true,
               python::arg("isomeric") = false, python::arg("kekulize") = false,
               python::arg("min_radius") = 1, python::arg("length") = 2048),
              "Creates a SECFP binary vector from a list of RDKit Mol "
              "instances."))
      .def("Distance", &MHFPEncoder::Distance, python::args("a", "b"))
      .staticmethod("Distance");
}

}  // namespace MHFPWrapper
}  // namespace RDKit