Package rdkit :: Package Chem :: Module rdForceFieldHelpers
[hide private]
[frames] | no frames]

Module rdForceFieldHelpers

Module containing functions to handle force fields

Functions [hide private]
 
GetUFFAngleBendParams(...)
GetUFFAngleBendParams( (Mol)mol, (int)idx1, (int)idx2, (int)idx3) -> object : Retrieves UFF angle bend parameters for atoms with indexes idx1, idx2, idx3 as a (ka, theta0) tuple, or None if no parameters could be found
 
GetUFFBondStretchParams(...)
GetUFFBondStretchParams( (Mol)mol, (int)idx1, (int)idx2) -> object : Retrieves UFF bond stretch parameters for atoms with indexes idx1, idx2 as a (kb, r0) tuple, or None if no parameters could be found
 
GetUFFInversionParams(...)
GetUFFInversionParams( (Mol)mol, (int)idx1, (int)idx2, (int)idx3, (int)idx4) -> object : Retrieves UFF inversion parameters for atoms with indexes idx1, idx2, idx3, idx4 as a K float value, or None if no parameters could be found
 
GetUFFTorsionParams(...)
GetUFFTorsionParams( (Mol)mol, (int)idx1, (int)idx2, (int)idx3, (int)idx4) -> object : Retrieves UFF torsion parameters for atoms with indexes idx1, idx2, idx3, idx4 as a V float value, or None if no parameters could be found
 
GetUFFVdWParams(...)
GetUFFVdWParams( (Mol)mol, (int)idx1, (int)idx2) -> object : Retrieves UFF van der Waals parameters for atoms with indexes idx1, idx2 as a (x_ij, D_ij) tuple, or None if no parameters could be found
 
MMFFGetMoleculeForceField(...)
MMFFGetMoleculeForceField( (Mol)mol, (MMFFMolProperties)pyMMFFMolProperties [, (float)nonBondedThresh=100.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]) -> ForceField : returns a MMFF force field for a molecule
 
MMFFGetMoleculeProperties(...)
MMFFGetMoleculeProperties( (Mol)mol [, (str)mmffVariant='MMFF94' [, (int)mmffVerbosity=0]]) -> MMFFMolProperties : returns a PyMMFFMolProperties object for a molecule, which is required by MMFFGetMoleculeForceField() and can be used to get/set MMFF properties
 
MMFFHasAllMoleculeParams(...)
MMFFHasAllMoleculeParams( (Mol)mol) -> bool : checks if MMFF parameters are available for all of a molecule's atoms
 
MMFFOptimizeMolecule(...)
MMFFOptimizeMolecule( (Mol)mol [, (str)mmffVariant='MMFF94' [, (int)maxIters=200 [, (float)nonBondedThresh=100.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]]) -> int : uses MMFF to optimize a molecule's structure
 
MMFFOptimizeMoleculeConfs(...)
MMFFOptimizeMoleculeConfs( (Mol)self [, (int)numThreads=1 [, (int)maxIters=200 [, (str)mmffVariant='MMFF94' [, (float)nonBondedThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]]]) -> object : uses MMFF to optimize all of a molecule's conformations
 
MMFFSanitizeMolecule(...)
MMFFSanitizeMolecule( (Mol)mol) -> int : sanitizes a molecule according to MMFF requirements.
 
UFFGetMoleculeForceField(...)
UFFGetMoleculeForceField( (Mol)mol [, (float)vdwThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]) -> ForceField : returns a UFF force field for a molecule
 
UFFHasAllMoleculeParams(...)
UFFHasAllMoleculeParams( (Mol)mol) -> bool : checks if UFF parameters are available for all of a molecule's atoms
 
UFFOptimizeMolecule(...)
UFFOptimizeMolecule( (Mol)self [, (int)maxIters=200 [, (float)vdwThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]) -> int : uses UFF to optimize a molecule's structure
 
UFFOptimizeMoleculeConfs(...)
UFFOptimizeMoleculeConfs( (Mol)self [, (int)numThreads=1 [, (int)maxIters=200 [, (float)vdwThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]]) -> object : uses UFF to optimize all of a molecule's conformations
Variables [hide private]
  __package__ = None
hash(x)
Function Details [hide private]

GetUFFAngleBendParams(...)

 

GetUFFAngleBendParams( (Mol)mol, (int)idx1, (int)idx2, (int)idx3) -> object :
    Retrieves UFF angle bend parameters for atoms with indexes idx1, idx2, idx3 as a (ka, theta0) tuple, or None if no parameters could be found

    C++ signature :
        _object* GetUFFAngleBendParams(RDKit::ROMol,unsigned int,unsigned int,unsigned int)

GetUFFBondStretchParams(...)

 

GetUFFBondStretchParams( (Mol)mol, (int)idx1, (int)idx2) -> object :
    Retrieves UFF bond stretch parameters for atoms with indexes idx1, idx2 as a (kb, r0) tuple, or None if no parameters could be found

    C++ signature :
        _object* GetUFFBondStretchParams(RDKit::ROMol,unsigned int,unsigned int)

GetUFFInversionParams(...)

 

GetUFFInversionParams( (Mol)mol, (int)idx1, (int)idx2, (int)idx3, (int)idx4) -> object :
    Retrieves UFF inversion parameters for atoms with indexes idx1, idx2, idx3, idx4 as a K float value, or None if no parameters could be found

    C++ signature :
        _object* GetUFFInversionParams(RDKit::ROMol,unsigned int,unsigned int,unsigned int,unsigned int)

GetUFFTorsionParams(...)

 

GetUFFTorsionParams( (Mol)mol, (int)idx1, (int)idx2, (int)idx3, (int)idx4) -> object :
    Retrieves UFF torsion parameters for atoms with indexes idx1, idx2, idx3, idx4 as a V float value, or None if no parameters could be found

    C++ signature :
        _object* GetUFFTorsionParams(RDKit::ROMol,unsigned int,unsigned int,unsigned int,unsigned int)

GetUFFVdWParams(...)

 

GetUFFVdWParams( (Mol)mol, (int)idx1, (int)idx2) -> object :
    Retrieves UFF van der Waals parameters for atoms with indexes idx1, idx2 as a (x_ij, D_ij) tuple, or None if no parameters could be found

    C++ signature :
        _object* GetUFFVdWParams(RDKit::ROMol,unsigned int,unsigned int)

MMFFGetMoleculeForceField(...)

 

MMFFGetMoleculeForceField( (Mol)mol, (MMFFMolProperties)pyMMFFMolProperties [, (float)nonBondedThresh=100.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]) -> ForceField :
    returns a MMFF force field for a molecule
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - pyMMFFMolProperties : PyMMFFMolProperties object as returned
                      by MMFFGetMoleculeProperties()
        - nonBondedThresh : used to exclude long-range non-bonded
                      interactions (defaults to 100.0)
        - confId : indicates which conformer to optimize
        - ignoreInterfragInteractions : if true, nonbonded terms between
                      fragments will not be added to the forcefield
    
    

    C++ signature :
        ForceFields::PyForceField* MMFFGetMoleculeForceField(RDKit::ROMol {lvalue},ForceFields::PyMMFFMolProperties* [,double=100.0 [,int=-1 [,bool=True]]])

MMFFGetMoleculeProperties(...)

 

MMFFGetMoleculeProperties( (Mol)mol [, (str)mmffVariant='MMFF94' [, (int)mmffVerbosity=0]]) -> MMFFMolProperties :
    returns a PyMMFFMolProperties object for a
      molecule, which is required by MMFFGetMoleculeForceField()
      and can be used to get/set MMFF properties
    
      
      ARGUMENTS:
    
        - mol : the molecule of interest
        - mmffVariant : "MMFF94" or "MMFF94s"
                      (defaults to "MMFF94")
        - mmffVerbosity : 0: none; 1: low; 2: high (defaults to 0).
    
    

    C++ signature :
        ForceFields::PyMMFFMolProperties* MMFFGetMoleculeProperties(RDKit::ROMol {lvalue} [,std::string='MMFF94' [,unsigned int=0]])

MMFFHasAllMoleculeParams(...)

 

MMFFHasAllMoleculeParams( (Mol)mol) -> bool :
    checks if MMFF parameters are available for all of a molecule's atoms
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
    
    

    C++ signature :
        bool MMFFHasAllMoleculeParams(RDKit::ROMol {lvalue})

MMFFOptimizeMolecule(...)

 

MMFFOptimizeMolecule( (Mol)mol [, (str)mmffVariant='MMFF94' [, (int)maxIters=200 [, (float)nonBondedThresh=100.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]]) -> int :
    uses MMFF to optimize a molecule's structure
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - mmffVariant : "MMFF94" or "MMFF94s"
        - maxIters : the maximum number of iterations (defaults to 200)
        - nonBondedThresh : used to exclude long-range non-bonded
                     interactions (defaults to 100.0)
        - confId : indicates which conformer to optimize
        - ignoreInterfragInteractions : if true, nonbonded terms between
                     fragments will not be added to the forcefield
    
     RETURNS: 0 if the optimization converged, -1 if the forcefield could
              not be set up, 1 if more iterations are required.
    
    

    C++ signature :
        int MMFFOptimizeMolecule(RDKit::ROMol {lvalue} [,std::string='MMFF94' [,int=200 [,double=100.0 [,int=-1 [,bool=True]]]]])

MMFFOptimizeMoleculeConfs(...)

 

MMFFOptimizeMoleculeConfs( (Mol)self [, (int)numThreads=1 [, (int)maxIters=200 [, (str)mmffVariant='MMFF94' [, (float)nonBondedThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]]]) -> object :
    uses MMFF to optimize all of a molecule's conformations
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - numThreads : the number of threads to use, only has an effect if the RDKit
                       was built with thread support (defaults to 1)
        - maxIters : the maximum number of iterations (defaults to 200)
        - mmffVariant : "MMFF94" or "MMFF94s"
        - nonBondedThresh : used to exclude long-range non-bonded
                      interactions (defaults to 100.0)
        - confId : indicates which conformer to optimize
        - ignoreInterfragInteractions : if true, nonbonded terms between
                      fragments will not be added to the forcefield.
    
     RETURNS: 0 if the optimization converged, 1 if more iterations are required.
    
    

    C++ signature :
        boost::python::api::object MMFFOptimizeMoleculeConfs(RDKit::ROMol {lvalue} [,unsigned int=1 [,int=200 [,std::string='MMFF94' [,double=10.0 [,int=-1 [,bool=True]]]]]])

MMFFSanitizeMolecule(...)

 

MMFFSanitizeMolecule( (Mol)mol) -> int :
    sanitizes a molecule according to MMFF requirements.
    
        - mol : the molecule of interest.
    
    

    C++ signature :
        unsigned int MMFFSanitizeMolecule(RDKit::ROMol {lvalue})

UFFGetMoleculeForceField(...)

 

UFFGetMoleculeForceField( (Mol)mol [, (float)vdwThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]) -> ForceField :
    returns a UFF force field for a molecule
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - vdwThresh : used to exclude long-range van der Waals interactions
                      (defaults to 10.0)
        - confId : indicates which conformer to optimize
        - ignoreInterfragInteractions : if true, nonbonded terms between
                      fragments will not be added to the forcefield.
    
    

    C++ signature :
        ForceFields::PyForceField* UFFGetMoleculeForceField(RDKit::ROMol {lvalue} [,double=10.0 [,int=-1 [,bool=True]]])

UFFHasAllMoleculeParams(...)

 

UFFHasAllMoleculeParams( (Mol)mol) -> bool :
    checks if UFF parameters are available for all of a molecule's atoms
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest.
    
    

    C++ signature :
        bool UFFHasAllMoleculeParams(RDKit::ROMol)

UFFOptimizeMolecule(...)

 

UFFOptimizeMolecule( (Mol)self [, (int)maxIters=200 [, (float)vdwThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]) -> int :
    uses UFF to optimize a molecule's structure
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - maxIters : the maximum number of iterations (defaults to 200)
        - vdwThresh : used to exclude long-range van der Waals interactions
                      (defaults to 10.0)
        - confId : indicates which conformer to optimize
        - ignoreInterfragInteractions : if true, nonbonded terms between
                      fragments will not be added to the forcefield.
    
     RETURNS: 0 if the optimization converged, 1 if more iterations are required.
    
    

    C++ signature :
        int UFFOptimizeMolecule(RDKit::ROMol {lvalue} [,int=200 [,double=10.0 [,int=-1 [,bool=True]]]])

UFFOptimizeMoleculeConfs(...)

 

UFFOptimizeMoleculeConfs( (Mol)self [, (int)numThreads=1 [, (int)maxIters=200 [, (float)vdwThresh=10.0 [, (int)confId=-1 [, (bool)ignoreInterfragInteractions=True]]]]]) -> object :
    uses UFF to optimize all of a molecule's conformations
    
     
     ARGUMENTS:
    
        - mol : the molecule of interest
        - numThreads : the number of threads to use, only has an effect if the RDKit
                       was built with thread support (defaults to 1)
        - maxIters : the maximum number of iterations (defaults to 200)
        - vdwThresh : used to exclude long-range van der Waals interactions
                      (defaults to 10.0)
        - confId : indicates which conformer to optimize
        - ignoreInterfragInteractions : if true, nonbonded terms between
                      fragments will not be added to the forcefield.
    
     RETURNS: 0 if the optimization converged, 1 if more iterations are required.
    
    

    C++ signature :
        boost::python::api::object UFFOptimizeMoleculeConfs(RDKit::ROMol {lvalue} [,unsigned int=1 [,int=200 [,double=10.0 [,int=-1 [,bool=True]]]]])