RDKit
Open-source cheminformatics and machine learning.
Loading...
Searching...
No Matches
RDDepict Namespace Reference

Classes

struct  Compute2DCoordParameters
 
struct  ConstrainedDepictionParams
 
class  CoordinateTemplates
 
class  DepictException
 
class  EmbeddedAtom
 Class that contains the data for an atoms that has already been embedded. More...
 
class  EmbeddedFrag
 Class containing a fragment of a molecule that has already been embedded. More...
 
struct  gtIIPair
 

Typedefs

typedef std::vector< const RDGeom::Point2D * > VECT_C_POINT
 
typedef std::pair< int, int > PAIR_I_I
 
typedef std::vector< PAIR_I_IVECT_PII
 
typedef std::priority_queue< PAIR_I_I, VECT_PII, gtIIPairPR_QUEUE
 
typedef std::pair< double, PAIR_I_IPAIR_D_I_I
 
typedef std::list< PAIR_D_I_ILIST_PAIR_DII
 
typedef std::pair< int, int > INT_PAIR
 
typedef std::vector< INT_PAIRINT_PAIR_VECT
 
typedef INT_PAIR_VECT::const_iterator INT_PAIR_VECT_CI
 
typedef std::pair< double, INT_PAIRDOUBLE_INT_PAIR
 
typedef boost::shared_array< double > DOUBLE_SMART_PTR
 
typedef std::map< unsigned int, EmbeddedAtomINT_EATOM_MAP
 
typedef INT_EATOM_MAP::iterator INT_EATOM_MAP_I
 
typedef INT_EATOM_MAP::const_iterator INT_EATOM_MAP_CI
 

Functions

RDKIT_CHEMREACTIONS_EXPORT void compute2DCoordsForReaction (RDKit::ChemicalReaction &rxn, double spacing=1.0, bool updateProps=true, bool canonOrient=false, unsigned int nFlipsPerSample=0, unsigned int nSamples=0, int sampleSeed=0, bool permuteDeg4Nodes=false)
 Generate 2D coordinates (a depiction) for a reaction.
 
RDKIT_DEPICTOR_EXPORT RDGeom::INT_POINT2D_MAP embedRing (const RDKit::INT_VECT &ring)
 Some utility functions used in generating 2D coordinates.
 
RDKIT_DEPICTOR_EXPORT void transformPoints (RDGeom::INT_POINT2D_MAP &nringCor, const RDGeom::Transform2D &trans)
 
RDKIT_DEPICTOR_EXPORT RDGeom::Point2D computeBisectPoint (const RDGeom::Point2D &rcr, double ang, const RDGeom::Point2D &nb1, const RDGeom::Point2D &nb2)
 Find a point that bisects the angle at rcr.
 
RDKIT_DEPICTOR_EXPORT void reflectPoints (RDGeom::INT_POINT2D_MAP &coordMap, const RDGeom::Point2D &loc1, const RDGeom::Point2D &loc2)
 Reflect a set of point through the line joining two point.
 
RDKIT_DEPICTOR_EXPORT RDGeom::Point2D reflectPoint (const RDGeom::Point2D &point, const RDGeom::Point2D &loc1, const RDGeom::Point2D &loc2)
 
RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT setNbrOrder (unsigned int aid, const RDKit::INT_VECT &nbrs, const RDKit::ROMol &mol)
 
RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT findNextRingToEmbed (const RDKit::INT_VECT &doneRings, const RDKit::VECT_INT_VECT &fusedRings, int &nextId)
 From a given set of rings find the ring the largest common elements with other rings.
 
template<class T >
RDKIT_DEPICTOR_EXPORTrankAtomsByRank (const RDKit::ROMol &mol, const T &commAtms, bool ascending=true)
 Sort a list of atoms by their CIP rank.
 
double computeSubAngle (unsigned int degree, RDKit::Atom::HybridizationType htype)
 computes a subangle for an atom of given hybridization and degree
 
int rotationDir (const RDGeom::Point2D &center, const RDGeom::Point2D &loc1, const RDGeom::Point2D &loc2, double remAngle)
 computes the rotation direction between two vectors
 
RDGeom::Point2D computeNormal (const RDGeom::Point2D &center, const RDGeom::Point2D &other)
 computes and return the normal of a vector between two points
 
double computeAngle (const RDGeom::Point2D &center, const RDGeom::Point2D &loc1, const RDGeom::Point2D &loc2)
 computes the rotation angle between two vectors
 
RDKIT_DEPICTOR_EXPORT int pickFirstRingToEmbed (const RDKit::ROMol &mol, const RDKit::VECT_INT_VECT &fusedRings)
 pick the ring to embed first in a fused system
 
RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT getRotatableBonds (const RDKit::ROMol &mol, unsigned int aid1, unsigned int aid2)
 find the rotatable bonds on the shortest path between two atoms we will ignore ring atoms, and double bonds which are marked cis/trans
 
RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT getAllRotatableBonds (const RDKit::ROMol &mol)
 find all the rotatable bonds in a molecule we will ignore ring atoms, and double bonds which are marked cis/trans
 
RDKIT_DEPICTOR_EXPORT void getNbrAtomAndBondIds (unsigned int aid, const RDKit::ROMol *mol, RDKit::INT_VECT &aids, RDKit::INT_VECT &bids)
 Get the ids of the atoms and bonds that are connected to aid.
 
RDKIT_DEPICTOR_EXPORT INT_PAIR_VECT findBondsPairsToPermuteDeg4 (const RDGeom::Point2D &center, const RDKit::INT_VECT &nbrBids, const VECT_C_POINT &nbrLocs)
 Find pairs of bonds that can be permuted at a non-ring degree 4 atom.
 
int getAtomDepictRank (const RDKit::Atom *at)
 returns the rank of the atom for determining draw order
 
RDKIT_DEPICTOR_EXPORT bool hasTerminalRGroupOrQueryHydrogen (const RDKit::ROMol &query)
 
RDKIT_DEPICTOR_EXPORT std::unique_ptr< RDKit::RWMolprepareTemplateForRGroups (RDKit::RWMol &templateMol)
 
RDKIT_DEPICTOR_EXPORT void reducedToFullMatches (const RDKit::RWMol &reducedQuery, const RDKit::RWMol &molHs, std::vector< RDKit::MatchVectType > &matches)
 
RDKIT_DEPICTOR_EXPORT bool invertWedgingIfMolHasFlipped (RDKit::ROMol &mol, const RDGeom::Transform3D &trans)
 
void RDKIT_DEPICTOR_EXPORT setRingSystemTemplates (const std::string templatePath)
 Set the path to the file containing the ring system templates.
 
void RDKIT_DEPICTOR_EXPORT addRingSystemTemplates (const std::string templatePath)
 Add ring system templates to be used in 2D coordinater generation. If there are duplicates, the most recently added template will be used.
 
void RDKIT_DEPICTOR_EXPORT loadDefaultRingSystemTemplates ()
 Load default ring system templates to be used in 2D coordinate generation.
 
RDKIT_DEPICTOR_EXPORT unsigned int compute2DCoords (RDKit::ROMol &mol, const Compute2DCoordParameters &params)
 Generate 2D coordinates (a depiction) for a molecule.
 
RDKIT_DEPICTOR_EXPORT unsigned int compute2DCoords (RDKit::ROMol &mol, const RDGeom::INT_POINT2D_MAP *coordMap=nullptr, bool canonOrient=false, bool clearConfs=true, unsigned int nFlipsPerSample=0, unsigned int nSamples=0, int sampleSeed=0, bool permuteDeg4Nodes=false, bool forceRDKit=false, bool useRingTemplates=false)
 Generate 2D coordinates (a depiction) for a molecule.
 
RDKIT_DEPICTOR_EXPORT unsigned int compute2DCoordsMimicDistMat (RDKit::ROMol &mol, const DOUBLE_SMART_PTR *dmat=nullptr, bool canonOrient=true, bool clearConfs=true, double weightDistMat=0.5, unsigned int nFlipsPerSample=3, unsigned int nSamples=100, int sampleSeed=25, bool permuteDeg4Nodes=true, bool forceRDKit=false)
 Compute the 2D coordinates such the interatom distances mimic those in a distance matrix.
 
RDKIT_DEPICTOR_EXPORT void generateDepictionMatching2DStructure (RDKit::ROMol &mol, const RDKit::ROMol &reference, const RDKit::MatchVectType &refMatchVect, int confId=-1, const ConstrainedDepictionParams &params=ConstrainedDepictionParams())
 Compute 2D coordinates where a piece of the molecule is constrained to have the same coordinates as a reference. Correspondences between reference and molecule atom indices are determined by refMatchVect.
 
RDKIT_DEPICTOR_EXPORT void generateDepictionMatching2DStructure (RDKit::ROMol &mol, const RDKit::ROMol &reference, const RDKit::MatchVectType &refMatchVect, int confId, bool forceRDKit)
 Overload.
 
RDKIT_DEPICTOR_EXPORT RDKit::MatchVectType generateDepictionMatching2DStructure (RDKit::ROMol &mol, const RDKit::ROMol &reference, int confId=-1, const RDKit::ROMol *referencePattern=static_cast< const RDKit::ROMol * >(nullptr), const ConstrainedDepictionParams &params=ConstrainedDepictionParams())
 Compute 2D coordinates constrained to a reference; the constraint can be hard (default) or soft.
 
RDKIT_DEPICTOR_EXPORT RDKit::MatchVectType generateDepictionMatching2DStructure (RDKit::ROMol &mol, const RDKit::ROMol &reference, int confId, const RDKit::ROMol *referencePattern, bool acceptFailure, bool forceRDKit=false, bool allowOptionalAttachments=false)
 Compute 2D coordinates where a piece of the molecule is constrained to have the same coordinates as a reference.
 
RDKIT_DEPICTOR_EXPORT void generateDepictionMatching3DStructure (RDKit::ROMol &mol, const RDKit::ROMol &reference, int confId=-1, RDKit::ROMol *referencePattern=nullptr, bool acceptFailure=false, bool forceRDKit=false)
 Generate a 2D depiction for a molecule where all or part of it mimics the coordinates of a 3D reference structure.
 
RDKIT_DEPICTOR_EXPORT void straightenDepiction (RDKit::ROMol &mol, int confId=-1, bool minimizeRotation=false)
 Rotate the 2D depiction such that the majority of bonds have an angle with the X axis which is a multiple of 30 degrees.
 
RDKIT_DEPICTOR_EXPORT double normalizeDepiction (RDKit::ROMol &mol, int confId=-1, int canonicalize=1, double scaleFactor=-1.0)
 Normalizes the 2D depiction.
 

Variables

RDKIT_DEPICTOR_EXPORT double BOND_LEN
 
RDKIT_DEPICTOR_EXPORT double COLLISION_THRES
 
RDKIT_DEPICTOR_EXPORT double BOND_THRES
 
RDKIT_DEPICTOR_EXPORT double ANGLE_OPEN
 
RDKIT_DEPICTOR_EXPORT unsigned int MAX_COLL_ITERS
 
RDKIT_DEPICTOR_EXPORT double HETEROATOM_COLL_SCALE
 
RDKIT_DEPICTOR_EXPORT unsigned int NUM_BONDS_FLIPS
 
RDKIT_DEPICTOR_EXPORT bool preferCoordGen
 

Typedef Documentation

◆ DOUBLE_INT_PAIR

typedef std::pair<double, INT_PAIR> RDDepict::DOUBLE_INT_PAIR

Definition at line 165 of file DepictUtils.h.

◆ DOUBLE_SMART_PTR

typedef boost::shared_array< double > RDDepict::DOUBLE_SMART_PTR

Definition at line 26 of file EmbeddedFrag.h.

◆ INT_EATOM_MAP

typedef std::map<unsigned int, EmbeddedAtom> RDDepict::INT_EATOM_MAP

Definition at line 132 of file EmbeddedFrag.h.

◆ INT_EATOM_MAP_CI

typedef INT_EATOM_MAP::const_iterator RDDepict::INT_EATOM_MAP_CI

Definition at line 134 of file EmbeddedFrag.h.

◆ INT_EATOM_MAP_I

typedef INT_EATOM_MAP::iterator RDDepict::INT_EATOM_MAP_I

Definition at line 133 of file EmbeddedFrag.h.

◆ INT_PAIR

typedef std::pair<int, int> RDDepict::INT_PAIR

Definition at line 161 of file DepictUtils.h.

◆ INT_PAIR_VECT

typedef std::vector<INT_PAIR> RDDepict::INT_PAIR_VECT

Definition at line 162 of file DepictUtils.h.

◆ INT_PAIR_VECT_CI

typedef INT_PAIR_VECT::const_iterator RDDepict::INT_PAIR_VECT_CI

Definition at line 163 of file DepictUtils.h.

◆ LIST_PAIR_DII

typedef std::list<PAIR_D_I_I> RDDepict::LIST_PAIR_DII

Definition at line 48 of file DepictUtils.h.

◆ PAIR_D_I_I

typedef std::pair<double, PAIR_I_I> RDDepict::PAIR_D_I_I

Definition at line 47 of file DepictUtils.h.

◆ PAIR_I_I

typedef std::pair<int, int> RDDepict::PAIR_I_I

Definition at line 37 of file DepictUtils.h.

◆ PR_QUEUE

typedef std::priority_queue<PAIR_I_I, VECT_PII, gtIIPair> RDDepict::PR_QUEUE

Definition at line 45 of file DepictUtils.h.

◆ VECT_C_POINT

typedef std::vector<const RDGeom::Point2D *> RDDepict::VECT_C_POINT

Definition at line 35 of file DepictUtils.h.

◆ VECT_PII

typedef std::vector<PAIR_I_I> RDDepict::VECT_PII

Definition at line 38 of file DepictUtils.h.

Function Documentation

◆ addRingSystemTemplates()

void RDKIT_DEPICTOR_EXPORT RDDepict::addRingSystemTemplates ( const std::string templatePath)

Add ring system templates to be used in 2D coordinater generation. If there are duplicates, the most recently added template will be used.

Parameters
templatePaththe file path to a file containing the ring system templates. Each template must be a single line in the file represented using CXSMILES, and the structure should be a single ring system.
Exceptions
DepictExceptionif any of the templates are invalid

◆ compute2DCoords() [1/2]

RDKIT_DEPICTOR_EXPORT unsigned int RDDepict::compute2DCoords ( RDKit::ROMol & mol,
const Compute2DCoordParameters & params )

Generate 2D coordinates (a depiction) for a molecule.

Parameters
molthe molecule were are interested in
paramsparameters used for 2D coordinate generation
Returns
ID of the conformation added to the molecule containing the 2D coordinates

References compute2DCoords().

Referenced by compute2DCoords(), and compute2DCoords().

◆ compute2DCoords() [2/2]

RDKIT_DEPICTOR_EXPORT unsigned int RDDepict::compute2DCoords ( RDKit::ROMol & mol,
const RDGeom::INT_POINT2D_MAP * coordMap = nullptr,
bool canonOrient = false,
bool clearConfs = true,
unsigned int nFlipsPerSample = 0,
unsigned int nSamples = 0,
int sampleSeed = 0,
bool permuteDeg4Nodes = false,
bool forceRDKit = false,
bool useRingTemplates = false )

Generate 2D coordinates (a depiction) for a molecule.

Parameters
molthe molecule were are interested in
coordMapa map of int to Point2D, between atom IDs and their locations. This is the container the user needs to fill if he/she wants to specify coordinates for a portion of the molecule, defaults to 0
canonOrientcanonicalize the orientation so that the long axes align with the x-axis etc.
clearConfsclear all existing conformations on the molecule before adding the 2D coordinates instead of simply adding to the list
nFlipsPerSample- the number of rotatable bonds that are flipped at random for each sample
nSamples- the number of samples
sampleSeed- seed for the random sampling process
permuteDeg4Nodes- try permuting the drawing order of bonds around atoms with four neighbors in order to improve the depiction
forceRDKit- use RDKit to generate coordinates even if preferCoordGen is set to true
useRingTemplateswhether to use ring system templates for generating initial coordinates
Returns
ID of the conformation added to the molecule containing the 2D coordinates

References compute2DCoords().

◆ compute2DCoordsForReaction()

RDKIT_CHEMREACTIONS_EXPORT void RDDepict::compute2DCoordsForReaction ( RDKit::ChemicalReaction & rxn,
double spacing = 1.0,
bool updateProps = true,
bool canonOrient = false,
unsigned int nFlipsPerSample = 0,
unsigned int nSamples = 0,
int sampleSeed = 0,
bool permuteDeg4Nodes = false )

Generate 2D coordinates (a depiction) for a reaction.

Parameters
rxnthe reaction we are interested in
spacingthe spacing between components of the reaction
updatePropsif set, properties such as conjugation and hybridization will be calculated for the reactant and product templates before generating coordinates. This should result in better depictions, but can lead to errors in some cases.
canonOrientcanonicalize the orientation so that the long axes align with the x-axis etc.
nFlipsPerSample- the number of rotatable bonds that are flipped at random for each sample
nSamples- the number of samples
sampleSeed- seed for the random sampling process
permuteDeg4Nodes- try permuting the drawing order of bonds around atoms with four neighbors in order to improve the depiction

for the other parameters see the documentation for compute2DCoords()

References compute2DCoordsForReaction().

Referenced by compute2DCoordsForReaction().

◆ compute2DCoordsMimicDistMat()

RDKIT_DEPICTOR_EXPORT unsigned int RDDepict::compute2DCoordsMimicDistMat ( RDKit::ROMol & mol,
const DOUBLE_SMART_PTR * dmat = nullptr,
bool canonOrient = true,
bool clearConfs = true,
double weightDistMat = 0.5,
unsigned int nFlipsPerSample = 3,
unsigned int nSamples = 100,
int sampleSeed = 25,
bool permuteDeg4Nodes = true,
bool forceRDKit = false )

Compute the 2D coordinates such the interatom distances mimic those in a distance matrix.

This function generates 2D coordinates such that the inter-atom distances mimic those specified via dmat. This is done by randomly sampling(flipping) the rotatable bonds in the molecule and evaluating a cost function which contains two components. The first component is the sum of inverse of the squared inter-atom distances, this helps in spreading the atoms far from each other. The second component is the sum of squares of the difference in distance between those in dmat and the generated structure. The user can adjust the relative importance of the two components via a adjustable parameter (see below)

ARGUMENTS:

Parameters
mol- molecule to generate coordinates for
dmat- the distance matrix we want to mimic, this is a symmetric N by N matrix where N is the number of atoms in mol. All negative entries in dmat are ignored.
canonOrient- canonicalize the orientation after the 2D embedding is done
clearConfs- clear any previously existing conformations on mol before adding a conformation
weightDistMat- A value between 0.0 and 1.0, this determines the importance of mimicing the inter atoms distances in dmat. (1.0 - weightDistMat) is the weight associated to spreading out the structure (density) in the cost function
nFlipsPerSample- the number of rotatable bonds that are flipped at random for each sample
nSamples- the number of samples
sampleSeed- seed for the random sampling process
permuteDeg4Nodes- try permuting the drawing order of bonds around atoms with four neighbors in order to improve the depiction
forceRDKit- use RDKit to generate coordinates even if preferCoordGen is set to true
Returns
ID of the conformation added to the molecule containing the 2D coordinates

References compute2DCoordsMimicDistMat().

Referenced by compute2DCoordsMimicDistMat().

◆ computeAngle()

double RDDepict::computeAngle ( const RDGeom::Point2D & center,
const RDGeom::Point2D & loc1,
const RDGeom::Point2D & loc2 )
inline

computes the rotation angle between two vectors

Parameters
centerthe common point
loc1endpoint 1
loc2endpoint 2
Returns
the angle (in radians)

Definition at line 273 of file DepictUtils.h.

References RDGeom::Point2D::angleTo().

◆ computeBisectPoint()

RDKIT_DEPICTOR_EXPORT RDGeom::Point2D RDDepict::computeBisectPoint ( const RDGeom::Point2D & rcr,
double ang,
const RDGeom::Point2D & nb1,
const RDGeom::Point2D & nb2 )

Find a point that bisects the angle at rcr.

The new point lies between nb1 and nb2. The line (rcr, newPt) bisects the angle 'ang' at rcr

◆ computeNormal()

RDGeom::Point2D RDDepict::computeNormal ( const RDGeom::Point2D & center,
const RDGeom::Point2D & other )
inline

computes and return the normal of a vector between two points

Parameters
centerthe common point
otherthe endpoint
Returns
the normal

Definition at line 256 of file DepictUtils.h.

References RDGeom::Point2D::normalize().

◆ computeSubAngle()

double RDDepict::computeSubAngle ( unsigned int degree,
RDKit::Atom::HybridizationType htype )
inline

computes a subangle for an atom of given hybridization and degree

Parameters
degreethe degree of the atom (number of neighbors)
htypethe atom's hybridization
Returns
the subangle (in radians)

Definition at line 185 of file DepictUtils.h.

References M_PI, RDKit::Atom::SP2, RDKit::Atom::SP3, and RDKit::Atom::UNSPECIFIED.

◆ embedRing()

RDKIT_DEPICTOR_EXPORT RDGeom::INT_POINT2D_MAP RDDepict::embedRing ( const RDKit::INT_VECT & ring)

Some utility functions used in generating 2D coordinates.

Embed a ring as a convex polygon in 2D

The process here is very straightforward:

We take the center of the ring to lie at the origin, so put the first point at the origin and then sweep anti-clockwise by an angle A = 360/n for the next point.

The length of the arm (l) we want to sweep is easy to compute given the bond length (b) we want to use for each bond in the ring (for now we will assume that this bond legnth is the same for all bonds in the ring:

l = b/sqrt(2*(1 - cos(A))

the above formula derives from the triangle formula, where side 'c' is given in terms of sides 'a' and 'b' as:

c = a^2 + b^2 - 2.a.b.cos(A)

where A is the angle between a and b

◆ findBondsPairsToPermuteDeg4()

RDKIT_DEPICTOR_EXPORT INT_PAIR_VECT RDDepict::findBondsPairsToPermuteDeg4 ( const RDGeom::Point2D & center,
const RDKit::INT_VECT & nbrBids,
const VECT_C_POINT & nbrLocs )

Find pairs of bonds that can be permuted at a non-ring degree 4 atom.

This function will return only those pairs that cannot be permuted by flipping a rotatble bond

   D
   |
   b3
   |

A-b1-B-b2-C | b4 | E For example in the above situation on the pairs (b1, b3) and (b1, b4) will be returned All other permutations can be achieved via a rotatable bond flip.

ARGUMENTS:

Parameters
center- location of the central atom
nbrBids- a vector (of length 4) containing the ids of the bonds to the neighbors
nbrLocs- locations of the neighbors

◆ findNextRingToEmbed()

RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT RDDepict::findNextRingToEmbed ( const RDKit::INT_VECT & doneRings,
const RDKit::VECT_INT_VECT & fusedRings,
int & nextId )

From a given set of rings find the ring the largest common elements with other rings.

◆ generateDepictionMatching2DStructure() [1/4]

RDKIT_DEPICTOR_EXPORT void RDDepict::generateDepictionMatching2DStructure ( RDKit::ROMol & mol,
const RDKit::ROMol & reference,
const RDKit::MatchVectType & refMatchVect,
int confId,
bool forceRDKit )

Overload.

ARGUMENTS:

Parameters
mol- the molecule to be aligned, this will come back with a single conformer.
reference- a molecule with the reference atoms to align to; this should have a depiction.
refMatchVect- a MatchVectType that will be used to generate the atom mapping between the molecule and the reference.
confId- the id of the reference conformation to use
forceRDKit- use RDKit to generate coordinates even if preferCoordGen is set to true

References generateDepictionMatching2DStructure().

◆ generateDepictionMatching2DStructure() [2/4]

RDKIT_DEPICTOR_EXPORT void RDDepict::generateDepictionMatching2DStructure ( RDKit::ROMol & mol,
const RDKit::ROMol & reference,
const RDKit::MatchVectType & refMatchVect,
int confId = -1,
const ConstrainedDepictionParams & params = ConstrainedDepictionParams() )

Compute 2D coordinates where a piece of the molecule is constrained to have the same coordinates as a reference. Correspondences between reference and molecule atom indices are determined by refMatchVect.

This function generates a depiction for a molecule where a piece of the molecule is constrained to have the same coordinates as a reference.

This is useful for, for example, generating depictions of SAR data sets so that the cores of the molecules are all oriented the same way. This overload allow to specify the (referenceAtom, molAtom) index pairs which should be matched as MatchVectType. Please note that the vector can be shorter than the number of atoms in the reference.

ARGUMENTS:

Parameters
mol- the molecule to be aligned, this will come back with a single conformer.
reference- a molecule with the reference atoms to align to; this should have a depiction.
refMatchVect- a MatchVectType that will be used to generate the atom mapping between the molecule and the reference.
confId- (optional) the id of the reference conformation to use
params- (optional) an instance of ConstrainedDepictionParams

References generateDepictionMatching2DStructure().

Referenced by generateDepictionMatching2DStructure(), generateDepictionMatching2DStructure(), generateDepictionMatching2DStructure(), and generateDepictionMatching2DStructure().

◆ generateDepictionMatching2DStructure() [3/4]

RDKIT_DEPICTOR_EXPORT RDKit::MatchVectType RDDepict::generateDepictionMatching2DStructure ( RDKit::ROMol & mol,
const RDKit::ROMol & reference,
int confId,
const RDKit::ROMol * referencePattern,
bool acceptFailure,
bool forceRDKit = false,
bool allowOptionalAttachments = false )

Compute 2D coordinates where a piece of the molecule is constrained to have the same coordinates as a reference.

This function generates a depiction for a molecule where a piece of the molecule is constrained to have the same coordinates as a reference.

This is useful, for example, for generating depictions of SAR data sets such that the cores of the molecules are all oriented the same way.

ARGUMENTS:

Parameters
mol- the molecule to be aligned, this will come back with a single conformer.
reference- a molecule with the reference atoms to align to; this should have a depiction.
confId- the id of the reference conformation to use
referencePattern- a query molecule to be used to generate the atom mapping between the molecule and the reference.
acceptFailure- if true, standard depictions will be generated for molecules that don't have a substructure match to the reference; if false, throws a DepictException.
forceRDKit- (optional) use RDKit to generate coordinates even if preferCoordGen is set to true
allowOptionalAttachments- (optional) if true, terminal dummy atoms in the reference are ignored if they match an implicit hydrogen in the molecule, and a constrained depiction is still attempted RETURNS:
Returns
MatchVectType with (queryAtomidx, molAtomIdx) pairs used for the constrained depiction

References generateDepictionMatching2DStructure().

◆ generateDepictionMatching2DStructure() [4/4]

RDKIT_DEPICTOR_EXPORT RDKit::MatchVectType RDDepict::generateDepictionMatching2DStructure ( RDKit::ROMol & mol,
const RDKit::ROMol & reference,
int confId = -1,
const RDKit::ROMol * referencePattern = static_cast< const RDKit::ROMol * >(nullptr),
const ConstrainedDepictionParams & params = ConstrainedDepictionParams() )

Compute 2D coordinates constrained to a reference; the constraint can be hard (default) or soft.

Hard (default, ConstrainedDepictionParams::alignOnly = false): Existing molecule coordinates, if present, are discarded; new coordinates are generated constraining a piece of the molecule to have the same coordinates as the reference, while the rest of the molecule is built around it. If ConstrainedDepictionParams::adjustMolBlockWedging is false (default), existing molblock wedging information is always preserved. If ConstrainedDepictionParams::adjustMolBlockWedging is true, existing molblock wedging information is preserved in case it only involves the invariant core and the core conformation has not changed, while it is cleared in case the wedging is also outside the invariant core, or core coordinates were changed. If ConstrainedDepictionParams::acceptFailure is set to true and no substructure match is found, coordinates will be recomputed from scratch, hence molblock wedging information will be cleared.

Soft (ConstrainedDepictionParams::alignOnly = true): Existing coordinates in the conformation identified by ConstrainedDepictionParams::existingConfId are preserved if present, otherwise unconstrained new coordinates are generated. Subsequently, coodinates undergo a rigid-body alignment to the reference. If ConstrainedDepictionParams::adjustMolBlockWedging is false (default), existing molblock wedging information is always preserved. If ConstrainedDepictionParams::adjustMolBlockWedging is true, existing molblock wedging information is inverted in case the rigid-body alignment involved a flip around the Z axis.

This is useful, for example, for generating depictions of SAR data sets such that the cores of the molecules are all oriented the same way.

ARGUMENTS:

Parameters
mol- the molecule to be aligned, this will come back with a single conformer.
reference- a molecule with the reference atoms to align to; this should have a depiction.
confId- (optional) the id of the reference conformation to use
referencePattern- (optional) a query molecule to be used to generate the atom mapping between the molecule and the reference.
params- (optional) a ConstrainedDepictionParams instance RETURNS:
Returns
MatchVectType with (queryAtomidx, molAtomIdx) pairs used for the constrained depiction

References generateDepictionMatching2DStructure().

◆ generateDepictionMatching3DStructure()

RDKIT_DEPICTOR_EXPORT void RDDepict::generateDepictionMatching3DStructure ( RDKit::ROMol & mol,
const RDKit::ROMol & reference,
int confId = -1,
RDKit::ROMol * referencePattern = nullptr,
bool acceptFailure = false,
bool forceRDKit = false )

Generate a 2D depiction for a molecule where all or part of it mimics the coordinates of a 3D reference structure.

Generates a depiction for a molecule where a piece of the molecule is constrained to have coordinates similar to those of a 3D reference structure.

ARGUMENTS:

Parameters
mol- the molecule to be aligned, this will come back with a single conformer containing 2D coordinates
reference- a molecule with the reference atoms to align to. By default this should be the same as mol, but with 3D coordinates
confId- (optional) the id of the reference conformation to use
refPattern- (optional) a query molecule to map a subset of the reference onto the mol, so that only some of the atoms are aligned.
acceptFailure- (optional) if true, standard depictions will be generated for molecules that don't match the reference or the referencePattern; if false, throws a DepictException.
forceRDKit- (optional) use RDKit to generate coordinates even if preferCoordGen is set to true

References generateDepictionMatching3DStructure().

Referenced by generateDepictionMatching3DStructure().

◆ getAllRotatableBonds()

RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT RDDepict::getAllRotatableBonds ( const RDKit::ROMol & mol)

find all the rotatable bonds in a molecule we will ignore ring atoms, and double bonds which are marked cis/trans

Note that rotatable in this context doesn't connect to the standard chemical definition of a rotatable bond; we're just talking about bonds than can be flipped in order to clean up the depiction.

Parameters
molthe molecule of interest
Returns
a set of the indices of the rotatable bonds

◆ getAtomDepictRank()

int RDDepict::getAtomDepictRank ( const RDKit::Atom * at)
inline

returns the rank of the atom for determining draw order

Definition at line 357 of file DepictUtils.h.

References RDKit::Atom::getAtomicNum(), and RDKit::Atom::getDegree().

◆ getNbrAtomAndBondIds()

RDKIT_DEPICTOR_EXPORT void RDDepict::getNbrAtomAndBondIds ( unsigned int aid,
const RDKit::ROMol * mol,
RDKit::INT_VECT & aids,
RDKit::INT_VECT & bids )

Get the ids of the atoms and bonds that are connected to aid.

◆ getRotatableBonds()

RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT RDDepict::getRotatableBonds ( const RDKit::ROMol & mol,
unsigned int aid1,
unsigned int aid2 )

find the rotatable bonds on the shortest path between two atoms we will ignore ring atoms, and double bonds which are marked cis/trans

Note that rotatable in this context doesn't connect to the standard chemical definition of a rotatable bond; we're just talking about bonds than can be flipped in order to clean up the depiction.

Parameters
molthe molecule of interest
aid1index of the first atom
aid2index of the second atom
Returns
a set of the indices of the rotatable bonds

◆ hasTerminalRGroupOrQueryHydrogen()

RDKIT_DEPICTOR_EXPORT bool RDDepict::hasTerminalRGroupOrQueryHydrogen ( const RDKit::ROMol & query)

◆ invertWedgingIfMolHasFlipped()

RDKIT_DEPICTOR_EXPORT bool RDDepict::invertWedgingIfMolHasFlipped ( RDKit::ROMol & mol,
const RDGeom::Transform3D & trans )

◆ loadDefaultRingSystemTemplates()

void RDKIT_DEPICTOR_EXPORT RDDepict::loadDefaultRingSystemTemplates ( )

Load default ring system templates to be used in 2D coordinate generation.

◆ normalizeDepiction()

RDKIT_DEPICTOR_EXPORT double RDDepict::normalizeDepiction ( RDKit::ROMol & mol,
int confId = -1,
int canonicalize = 1,
double scaleFactor = -1.0 )

Normalizes the 2D depiction.

If canonicalize is != 0, the depiction is subjected to a canonical transformation such that its main axis is aligned along the X axis (canonicalize >0, the default) or the Y axis (canonicalize <0). If canonicalize is 0, no canonicalization takes place. If scaleFactor is <0.0 (the default) the depiction is scaled such that bond lengths conform to RDKit standards. The applied scaling factor is returned.

ARGUMENTS:

Parameters
mol- the molecule to be normalized
confId- (optional) the id of the reference conformation to use
canonicalize- (optional) if != 0, a canonical transformation is applied: if >0 (the default), the main molecule axis is aligned to the X axis, if <0 to the Y axis. If 0, no canonical transformation is applied.
scaleFactor- (optional) if >0.0, the scaling factor to apply. The default (-1.0) means that the depiction is automatically scaled such that bond lengths are the standard RDKit ones. RETURNS:
Returns
the applied scaling factor.

References normalizeDepiction().

Referenced by normalizeDepiction().

◆ pickFirstRingToEmbed()

RDKIT_DEPICTOR_EXPORT int RDDepict::pickFirstRingToEmbed ( const RDKit::ROMol & mol,
const RDKit::VECT_INT_VECT & fusedRings )

pick the ring to embed first in a fused system

Parameters
molthe molecule of interest
fusedRingsthe collection of the molecule's fused rings
Returns
the index of the ring with the least number of substitutions

◆ prepareTemplateForRGroups()

RDKIT_DEPICTOR_EXPORT std::unique_ptr< RDKit::RWMol > RDDepict::prepareTemplateForRGroups ( RDKit::RWMol & templateMol)

◆ rankAtomsByRank()

template<class T >
RDKIT_DEPICTOR_EXPORT T RDDepict::rankAtomsByRank ( const RDKit::ROMol & mol,
const T & commAtms,
bool ascending = true )

Sort a list of atoms by their CIP rank.

Parameters
molmolecule of interest
commAtmsatoms that need to be ranked
ascendingsort to an ascending order or a descending order

◆ reducedToFullMatches()

RDKIT_DEPICTOR_EXPORT void RDDepict::reducedToFullMatches ( const RDKit::RWMol & reducedQuery,
const RDKit::RWMol & molHs,
std::vector< RDKit::MatchVectType > & matches )

◆ reflectPoint()

RDKIT_DEPICTOR_EXPORT RDGeom::Point2D RDDepict::reflectPoint ( const RDGeom::Point2D & point,
const RDGeom::Point2D & loc1,
const RDGeom::Point2D & loc2 )

◆ reflectPoints()

RDKIT_DEPICTOR_EXPORT void RDDepict::reflectPoints ( RDGeom::INT_POINT2D_MAP & coordMap,
const RDGeom::Point2D & loc1,
const RDGeom::Point2D & loc2 )

Reflect a set of point through the line joining two point.

ARGUMENTS:

Parameters
coordMapa map of <int, point2D> going from atom id to current coordinates of the points that need to be reflected: The coordinates are overwritten
loc1the first point of the line that is to be used as a mirror
loc2the second point of the line to be used as a mirror

◆ rotationDir()

int RDDepict::rotationDir ( const RDGeom::Point2D & center,
const RDGeom::Point2D & loc1,
const RDGeom::Point2D & loc2,
double remAngle )
inline

computes the rotation direction between two vectors

Let:

v1 = loc1 - center

v2 = loc2 - center

If remaining angle(v1, v2) is < 180 and corss(v1, v2) > 0.0 then the rotation dir is +1.0

else if remAngle(v1, v2) is > 180 and cross(v1, v2) < 0.0 then rotation dir is -1.0

else if remAngle(v1, v2) is < 180 and cross(v1, v2) < 0.0 then rotation dir is -1.0

finally if remAngle(v1, v2) is > 180 and cross(v1, v2) < 0.0 then rotation dir is +1.0

Parameters
centerthe common point
loc1endpoint 1
loc2endpoint 2
remAnglethe remaining angle about center in radians
Returns
the rotation direction (1 or -1)

Definition at line 234 of file DepictUtils.h.

References M_PI, and RDGeom::Point2D::x.

◆ setNbrOrder()

RDKIT_DEPICTOR_EXPORT RDKit::INT_VECT RDDepict::setNbrOrder ( unsigned int aid,
const RDKit::INT_VECT & nbrs,
const RDKit::ROMol & mol )

Set the neighbors yet to added to aid such that the atoms with the most subs fall on opposite sides

Ok this needs some explanation

  • Let A, B, C, D be the substituent on the central atom X (given by atom index aid)
  • also let A be the atom that is already embedded
  • Now we want the order in which the remaining neighbors B,C,D are added to X such that the atoms with atom with largest number of substituent fall on opposite sides of X so as to minimize atom clashes later in the depiction

E.g. let say we have the following situation

          B
       |  |
       A--X--C
       |  |
        --D--
          |

In this case the number substituent of A, B, C, D are 3, 1, 1, 4 respectively so want to A and D to go opposite sides and so that we draw

          B
       |  |  |
       A--X--D--
       |  |  |
          C

And the correct ordering of the neighbors is B,D,C

◆ setRingSystemTemplates()

void RDKIT_DEPICTOR_EXPORT RDDepict::setRingSystemTemplates ( const std::string templatePath)

Set the path to the file containing the ring system templates.

Parameters
templatePaththe file path to a file containing the ring system templates. Each template must be a single line in the file represented using CXSMILES, and the structure should be a single ring system.
Exceptions
DepictExceptionif any of the templates are invalid

◆ straightenDepiction()

RDKIT_DEPICTOR_EXPORT void RDDepict::straightenDepiction ( RDKit::ROMol & mol,
int confId = -1,
bool minimizeRotation = false )

Rotate the 2D depiction such that the majority of bonds have an angle with the X axis which is a multiple of 30 degrees.

ARGUMENTS:

Parameters
mol- the molecule to be rotated
confId- (optional) the id of the reference conformation to use
minimizeRotation- (optional) if false (the default), the molecule is rotated such that the majority of bonds have an angle with the X axis of 30 or 90 degrees. If true, the minimum rotation is applied such that the majority of bonds have an angle with the X axis of 0, 30, 60, or 90 degrees, with the goal of altering the initial orientation as little as possible .

References straightenDepiction().

Referenced by straightenDepiction().

◆ transformPoints()

RDKIT_DEPICTOR_EXPORT void RDDepict::transformPoints ( RDGeom::INT_POINT2D_MAP & nringCor,
const RDGeom::Transform2D & trans )

Variable Documentation

◆ ANGLE_OPEN

RDKIT_DEPICTOR_EXPORT double RDDepict::ANGLE_OPEN
extern

◆ BOND_LEN

RDKIT_DEPICTOR_EXPORT double RDDepict::BOND_LEN
extern

◆ BOND_THRES

RDKIT_DEPICTOR_EXPORT double RDDepict::BOND_THRES
extern

◆ COLLISION_THRES

RDKIT_DEPICTOR_EXPORT double RDDepict::COLLISION_THRES
extern

◆ HETEROATOM_COLL_SCALE

RDKIT_DEPICTOR_EXPORT double RDDepict::HETEROATOM_COLL_SCALE
extern

◆ MAX_COLL_ITERS

RDKIT_DEPICTOR_EXPORT unsigned int RDDepict::MAX_COLL_ITERS
extern

◆ NUM_BONDS_FLIPS

RDKIT_DEPICTOR_EXPORT unsigned int RDDepict::NUM_BONDS_FLIPS
extern

◆ preferCoordGen

RDKIT_DEPICTOR_EXPORT bool RDDepict::preferCoordGen
extern