Skip to content

Commit

Permalink
- remove deprecation warnings by switching MinimalLib to fingerprint …
Browse files Browse the repository at this point in the history
…generators (rdkit#7938)

- add missing Fingerprints dependency to MinimalLib CMakeLists.txt
- remove conditional RGroupDecomposition dependency in MinimalLib CMakeLists.txt, as it is always needed because of relabelMappedDummies()

wip

Co-authored-by: ptosco <[email protected]>
  • Loading branch information
ptosco and ptosco authored Oct 30, 2024
1 parent eacc365 commit 64061b6
Show file tree
Hide file tree
Showing 2 changed files with 49 additions and 24 deletions.
10 changes: 3 additions & 7 deletions Code/MinimalLib/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,7 @@ if(RDK_BUILD_MINIMAL_LIB)
"CIPLabeler_static;MolDraw2D_static;Depictor_static;"
"Descriptors_static;SubstructMatch_static;FileParsers_static;"
"SmilesParse_static;GraphMol_static;RDGeometryLib_static;"
"RDGeneral_static;RGroupDecomposition_static")
"RDGeneral_static;RGroupDecomposition_static;Fingerprints_static")
if(RDK_BUILD_INCHI_SUPPORT)
add_definitions(-DRDK_BUILD_INCHI_SUPPORT)
set(MINIMAL_LIB_LIBRARIES "${MINIMAL_LIB_LIBRARIES};RDInchiLib_static")
Expand All @@ -23,10 +23,6 @@ if(RDK_BUILD_MINIMAL_LIB)
add_definitions(-DRDK_BUILD_MINIMAL_LIB_MCS)
set(MINIMAL_LIB_LIBRARIES "${MINIMAL_LIB_LIBRARIES};FMCS_static")
endif()
if(RDK_BUILD_MINIMAL_LIB_RGROUPDECOMP)
add_definitions(-DRDK_BUILD_MINIMAL_LIB_RGROUPDECOMP)
set(MINIMAL_LIB_LIBRARIES "${MINIMAL_LIB_LIBRARIES};RGroupDecomposition")
endif()
if(RDK_BUILD_MINIMAL_LIB_MMPA)
add_definitions(-DRDK_BUILD_MINIMAL_LIB_MMPA)
set(MINIMAL_LIB_LIBRARIES "${MINIMAL_LIB_LIBRARIES};MMPA_static")
Expand All @@ -52,8 +48,8 @@ if(RDK_BUILD_CFFI_LIB)
ForceField Alignment
MolInterchange Abbreviations CIPLabeler
MolDraw2D Depictor Descriptors
SubstructMatch FileParsers
SmilesParse GraphMol RDGeometryLib RDGeneral RGroupDecomposition)
SubstructMatch FileParsers SmilesParse GraphMol
RDGeometryLib RDGeneral RGroupDecomposition Fingerprints)
if(RDK_BUILD_INCHI_SUPPORT)
add_definitions(-DRDK_BUILD_INCHI_SUPPORT)
list(APPEND LIBS_TO_USE RDInchiLib)
Expand Down
63 changes: 46 additions & 17 deletions Code/MinimalLib/common.h
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,9 @@
// of the RDKit source tree.
//
#pragma once
#ifdef __GNUC__
#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
#endif

#include <string>
#include <RDGeneral/versions.h>
Expand All @@ -26,7 +28,10 @@
#include <GraphMol/Descriptors/MolDescriptors.h>
#include <GraphMol/Fingerprints/Fingerprints.h>
#include <GraphMol/Fingerprints/MorganFingerprints.h>
#include <GraphMol/Fingerprints/AtomPairs.h>
#include <GraphMol/Fingerprints/AtomPairGenerator.h>
#include <GraphMol/Fingerprints/TopologicalTorsionGenerator.h>
#include <GraphMol/Fingerprints/MorganGenerator.h>
#include <GraphMol/Fingerprints/RDKitFPGenerator.h>
#include <GraphMol/Fingerprints/MACCS.h>
#ifdef RDK_BUILD_AVALON_SUPPORT
#include <External/AvalonTools/AvalonTools.h>
Expand Down Expand Up @@ -844,8 +849,9 @@ std::unique_ptr<RWMol> do_fragment_parent(RWMol &mol,

std::unique_ptr<ExplicitBitVect> morgan_fp_as_bitvect(
const RWMol &mol, const char *details_json) {
size_t radius = 2;
size_t nBits = 2048;
unsigned int radius = 2;
unsigned int nBits = 2048;
bool useCountSimulation = false;
bool useChirality = false;
bool useBondTypes = true;
bool includeRedundantEnvironments = false;
Expand All @@ -858,26 +864,31 @@ std::unique_ptr<ExplicitBitVect> morgan_fp_as_bitvect(
boost::property_tree::read_json(ss, pt);
LPT_OPT_GET(radius);
LPT_OPT_GET(nBits);
LPT_OPT_GET(useCountSimulation);
LPT_OPT_GET(useChirality);
LPT_OPT_GET(useBondTypes);
LPT_OPT_GET(includeRedundantEnvironments);
LPT_OPT_GET(onlyNonzeroInvariants);
}
auto fp = MorganFingerprints::getFingerprintAsBitVect(
mol, radius, nBits, nullptr, nullptr, useChirality, useBondTypes,
onlyNonzeroInvariants, nullptr, includeRedundantEnvironments);
return std::unique_ptr<ExplicitBitVect>{fp};
std::unique_ptr<FingerprintGenerator<std::uint32_t>> morganFPGen{
RDKit::MorganFingerprint::getMorganGenerator<std::uint32_t>(
radius, useCountSimulation, useChirality, useBondTypes,
onlyNonzeroInvariants, includeRedundantEnvironments, nullptr, nullptr,
nBits)};
return std::unique_ptr<ExplicitBitVect>(morganFPGen->getFingerprint(mol));
}

std::unique_ptr<ExplicitBitVect> rdkit_fp_as_bitvect(const RWMol &mol,
const char *details_json) {
static const std::vector<std::uint32_t> countBounds{1, 2, 4, 8};
unsigned int minPath = 1;
unsigned int maxPath = 7;
unsigned int nBits = 2048;
unsigned int nBitsPerHash = 2;
bool useHs = true;
bool branchedPaths = true;
bool useBondOrder = true;
bool useCountSimulation = false;
if (details_json && strlen(details_json)) {
// FIX: this should eventually be moved somewhere else
std::istringstream ss;
Expand All @@ -891,10 +902,13 @@ std::unique_ptr<ExplicitBitVect> rdkit_fp_as_bitvect(const RWMol &mol,
LPT_OPT_GET(useHs);
LPT_OPT_GET(branchedPaths);
LPT_OPT_GET(useBondOrder);
LPT_OPT_GET(useCountSimulation);
}
auto fp = RDKFingerprintMol(mol, minPath, maxPath, nBits, nBitsPerHash, useHs,
0, 128, branchedPaths, useBondOrder);
return std::unique_ptr<ExplicitBitVect>{fp};
std::unique_ptr<FingerprintGenerator<std::uint32_t>> rdkFPGen{
RDKit::RDKitFP::getRDKitFPGenerator<std::uint32_t>(
minPath, maxPath, useHs, branchedPaths, useBondOrder, nullptr,
useCountSimulation, countBounds, nBits, nBitsPerHash)};
return std::unique_ptr<ExplicitBitVect>(rdkFPGen->getFingerprint(mol));
}

std::unique_ptr<ExplicitBitVect> pattern_fp_as_bitvect(
Expand All @@ -918,24 +932,34 @@ std::unique_ptr<ExplicitBitVect> pattern_fp_as_bitvect(
std::unique_ptr<ExplicitBitVect> topological_torsion_fp_as_bitvect(
const RWMol &mol, const char *details_json) {
unsigned int nBits = 2048;
unsigned int torsionAtomCount = 4;
bool useChirality = false;
bool useCountSimulation = true;
if (details_json && strlen(details_json)) {
// FIX: this should eventually be moved somewhere else
std::istringstream ss;
ss.str(details_json);
boost::property_tree::ptree pt;
boost::property_tree::read_json(ss, pt);
LPT_OPT_GET(nBits);
LPT_OPT_GET(torsionAtomCount);
LPT_OPT_GET(useChirality);
LPT_OPT_GET(useCountSimulation);
}
auto fp =
AtomPairs::getHashedTopologicalTorsionFingerprintAsBitVect(mol, nBits);
return std::unique_ptr<ExplicitBitVect>{fp};
std::unique_ptr<FingerprintGenerator<std::uint64_t>> topoTorsFPGen{
RDKit::TopologicalTorsion::getTopologicalTorsionGenerator<std::uint64_t>(
useChirality, torsionAtomCount, nullptr, useCountSimulation, nBits)};
return std::unique_ptr<ExplicitBitVect>(topoTorsFPGen->getFingerprint(mol));
}

std::unique_ptr<ExplicitBitVect> atom_pair_fp_as_bitvect(
const RWMol &mol, const char *details_json) {
unsigned int nBits = 2048;
unsigned int minLength = 1;
unsigned int maxLength = 30;
bool useChirality = false;
bool use2D = true;
bool useCountSimulation = true;
if (details_json && strlen(details_json)) {
// FIX: this should eventually be moved somewhere else
std::istringstream ss;
Expand All @@ -945,10 +969,15 @@ std::unique_ptr<ExplicitBitVect> atom_pair_fp_as_bitvect(
LPT_OPT_GET(nBits);
LPT_OPT_GET(minLength);
LPT_OPT_GET(maxLength);
}
auto fp = AtomPairs::getHashedAtomPairFingerprintAsBitVect(
mol, nBits, minLength, maxLength);
return std::unique_ptr<ExplicitBitVect>{fp};
LPT_OPT_GET(useChirality);
LPT_OPT_GET(use2D);
LPT_OPT_GET(useCountSimulation);
}
std::unique_ptr<FingerprintGenerator<std::uint32_t>> atomPairFPGen{
RDKit::AtomPair::getAtomPairGenerator<std::uint32_t>(
minLength, maxLength, useChirality, use2D, nullptr,
useCountSimulation, nBits)};
return std::unique_ptr<ExplicitBitVect>(atomPairFPGen->getFingerprint(mol));
}

std::unique_ptr<ExplicitBitVect> maccs_fp_as_bitvect(const RWMol &mol) {
Expand Down

0 comments on commit 64061b6

Please sign in to comment.