get.fingerprint {rcdk} | R Documentation |
Generate molecular fingerprints
Description
'get.fingerprint' returns a 'fingerprint' object representing molecular fingerprint of the input molecule.
Usage
get.fingerprint(
molecule,
type = "standard",
fp.mode = "bit",
depth = 6,
size = 1024,
substructure.pattern = character(),
circular.type = "ECFP6",
verbose = FALSE
)
Arguments
molecule |
A |
type |
The type of fingerprint. Possible values are:
|
fp.mode |
The style of fingerprint. Specifying "'bit'" will return a binary fingerprint, "'raw'" returns the the original representation (usually sequence of integers) and "'count'" returns the fingerprint as a sequence of counts. |
depth |
The search depth. This argument is ignored for the 'pubchem', 'maccs', 'kr' and 'estate' fingerprints |
size |
The final length of the fingerprint. This argument is ignored for the 'pubchem', 'maccs', 'kr', 'signature', 'circular' and 'estate' fingerprints |
substructure.pattern |
List of characters containing the SMARTS pattern to match. If the an empty list is provided (default) than the functional groups substructures (default in CDK) are used. |
circular.type |
Name of the circular fingerprint type that should be computed given as string. Possible values are: 'ECFP0', 'ECFP2', 'ECFP4', 'ECFP6' (default), 'FCFP0', 'FCFP2', 'FCFP4' and 'FCFP6'. |
verbose |
Verbose output if |
Value
an S4 object of class fingerprint-class
or featvec-class
,
which can be manipulated with the fingerprint package.
Author(s)
Rajarshi Guha (rajarshi.guha@gmail.com)
Examples
## get some molecules
sp <- get.smiles.parser()
smiles <- c('CCC', 'CCN', 'CCN(C)(C)', 'c1ccccc1Cc1ccccc1','C1CCC1CC(CN(C)(C))CC(=O)CC')
mols <- parse.smiles(smiles)
## get a single fingerprint using the standard
## (hashed, path based) fingerprinter
fp <- get.fingerprint(mols[[1]])
## get MACCS keys for all the molecules
fps <- lapply(mols, get.fingerprint, type='maccs')
## get Signature fingerprint
## feature, count fingerprinter
fps <- lapply(mols, get.fingerprint, type='signature', fp.mode='raw')
## get Substructure fingerprint for functional group fragments
fps <- lapply(mols, get.fingerprint, type='substructure')
## get Substructure count fingerprint for user defined fragments
mol1 <- parse.smiles("c1ccccc1CCC")[[1]]
smarts <- c("c1ccccc1", "[CX4H3][#6]", "[CX2]#[CX2]")
fps <- get.fingerprint(mol1, type='substructure', fp.mode='count',
substructure.pattern=smarts)
## get ECFP0 count fingerprints
mol2 <- parse.smiles("C1=CC=CC(=C1)CCCC2=CC=CC=C2")[[1]]
fps <- get.fingerprint(mol2, type='circular', fp.mode='count', circular.type='ECFP0')